particle.C 32.2 KB
Newer Older
andy's avatar
andy committed
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2017 OpenFOAM Foundation
6
     \\/     M anipulation  | Copyright (C) 2018 OpenCFD Ltd.
andy's avatar
andy committed
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
-------------------------------------------------------------------------------
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/>.

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

#include "particle.H"
#include "transform.H"
28
29
#include "treeDataCell.H"
#include "cubicEqn.H"
30
#include "registerSwitch.H"
andy's avatar
andy committed
31
32
33

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

34
35
36
37
namespace Foam
{
    defineTypeNameAndDebug(particle, 0);
}
andy's avatar
andy committed
38

39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
const Foam::scalar Foam::particle::negativeSpaceDisplacementFactor = 1.01;

Foam::label Foam::particle::particleCount_ = 0;

bool Foam::particle::writeLagrangianCoordinates = true;

bool Foam::particle::writeLagrangianPositions
(
    Foam::debug::infoSwitch("writeLagrangianPositions", 0)
);

registerInfoSwitch
(
    "writeLagrangianPositions",
    bool,
    Foam::particle::writeLagrangianPositions
);

57

58
// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
59

60
void Foam::particle::stationaryTetGeometry
61
62
63
64
65
66
67
(
    vector& centre,
    vector& base,
    vector& vertex1,
    vector& vertex2
) const
{
68
69
70
71
72
73
74
75
    const triFace triIs(currentTetIndices().faceTriIs(mesh_));
    const vectorField& ccs = mesh_.cellCentres();
    const pointField& pts = mesh_.points();

    centre = ccs[celli_];
    base = pts[triIs[0]];
    vertex1 = pts[triIs[1]];
    vertex2 = pts[triIs[2]];
76
77
78
}


79
Foam::barycentricTensor Foam::particle::stationaryTetTransform() const
80
81
{
    vector centre, base, vertex1, vertex2;
82
    stationaryTetGeometry(centre, base, vertex1, vertex2);
83
84
85
86
87

    return barycentricTensor(centre, base, vertex1, vertex2);
}


88
void Foam::particle::stationaryTetReverseTransform
89
90
91
92
93
94
(
    vector& centre,
    scalar& detA,
    barycentricTensor& T
) const
{
95
    barycentricTensor A = stationaryTetTransform();
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

    const vector ab = A.b() - A.a();
    const vector ac = A.c() - A.a();
    const vector ad = A.d() - A.a();
    const vector bc = A.c() - A.b();
    const vector bd = A.d() - A.b();

    centre = A.a();

    detA = ab & (ac ^ ad);

    T = barycentricTensor
    (
        bd ^ bc,
        ac ^ ad,
        ad ^ ab,
        ab ^ ac
    );
}


void Foam::particle::movingTetGeometry
(
    const scalar fraction,
    Pair<vector>& centre,
    Pair<vector>& base,
    Pair<vector>& vertex1,
    Pair<vector>& vertex2
) const
{
126
    const triFace triIs(currentTetIndices().faceTriIs(mesh_));
127
128
129
130
131
132
133
134
135
136
137
138
139
    const pointField& ptsOld = mesh_.oldPoints();
    const pointField& ptsNew = mesh_.points();

    // !!! <-- We would be better off using mesh_.cellCentres() here. However,
    // we need to put a mesh_.oldCellCentres() method in for this to work. The
    // values obtained from the mesh and those obtained from the cell do not
    // necessarily match. See mantis #1993.
    const vector ccOld = mesh_.cells()[celli_].centre(ptsOld, mesh_.faces());
    const vector ccNew = mesh_.cells()[celli_].centre(ptsNew, mesh_.faces());

    // Old and new points and cell centres are not sub-cycled. If we are sub-
    // cycling, then we have to account for the timestep change here by
    // modifying the fractions that we take of the old and new geometry.
140
141
    const Pair<scalar> s = stepFractionSpan();
    const scalar f0 = s[0] + stepFraction_*s[1], f1 = fraction*s[1];
142
143

    centre[0] = ccOld + f0*(ccNew - ccOld);
144
145
146
    base[0] = ptsOld[triIs[0]] + f0*(ptsNew[triIs[0]] - ptsOld[triIs[0]]);
    vertex1[0] = ptsOld[triIs[1]] + f0*(ptsNew[triIs[1]] - ptsOld[triIs[1]]);
    vertex2[0] = ptsOld[triIs[2]] + f0*(ptsNew[triIs[2]] - ptsOld[triIs[2]]);
147
148

    centre[1] = f1*(ccNew - ccOld);
149
150
151
    base[1] = f1*(ptsNew[triIs[0]] - ptsOld[triIs[0]]);
    vertex1[1] = f1*(ptsNew[triIs[1]] - ptsOld[triIs[1]]);
    vertex2[1] = f1*(ptsNew[triIs[2]] - ptsOld[triIs[2]]);
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
}


Foam::Pair<Foam::barycentricTensor> Foam::particle::movingTetTransform
(
    const scalar fraction
) const
{
    Pair<vector> centre, base, vertex1, vertex2;
    movingTetGeometry(fraction, centre, base, vertex1, vertex2);

    return
        Pair<barycentricTensor>
        (
            barycentricTensor(centre[0], base[0], vertex1[0], vertex2[0]),
            barycentricTensor(centre[1], base[1], vertex1[1], vertex2[1])
        );
}


void Foam::particle::movingTetReverseTransform
(
    const scalar fraction,
    Pair<vector>& centre,
    FixedList<scalar, 4>& detA,
    FixedList<barycentricTensor, 3>& T
) const
{
    Pair<barycentricTensor> A = movingTetTransform(fraction);

    const Pair<vector> ab(A[0].b() - A[0].a(), A[1].b() - A[1].a());
    const Pair<vector> ac(A[0].c() - A[0].a(), A[1].c() - A[1].a());
    const Pair<vector> ad(A[0].d() - A[0].a(), A[1].d() - A[1].a());
    const Pair<vector> bc(A[0].c() - A[0].b(), A[1].c() - A[1].b());
    const Pair<vector> bd(A[0].d() - A[0].b(), A[1].d() - A[1].b());

    centre[0] = A[0].a();
    centre[1] = A[1].a();

    detA[0] = ab[0] & (ac[0] ^ ad[0]);
    detA[1] =
        (ab[1] & (ac[0] ^ ad[0]))
      + (ab[0] & (ac[1] ^ ad[0]))
      + (ab[0] & (ac[0] ^ ad[1]));
    detA[2] =
        (ab[0] & (ac[1] ^ ad[1]))
      + (ab[1] & (ac[0] ^ ad[1]))
      + (ab[1] & (ac[1] ^ ad[0]));
    detA[3] = ab[1] & (ac[1] ^ ad[1]);

    T[0] = barycentricTensor
    (
        bd[0] ^ bc[0],
        ac[0] ^ ad[0],
        ad[0] ^ ab[0],
        ab[0] ^ ac[0]
    );
    T[1] = barycentricTensor
    (
        (bd[0] ^ bc[1]) + (bd[1] ^ bc[0]),
        (ac[0] ^ ad[1]) + (ac[1] ^ ad[0]),
        (ad[0] ^ ab[1]) + (ad[1] ^ ab[0]),
        (ab[0] ^ ac[1]) + (ab[1] ^ ac[0])
    );
    T[2] = barycentricTensor
    (
        bd[1] ^ bc[1],
        ac[1] ^ ad[1],
        ad[1] ^ ab[1],
        ab[1] ^ ac[1]
    );
}


void Foam::particle::reflect()
{
    Swap(coordinates_.c(), coordinates_.d());
}


void Foam::particle::rotate(const bool reverse)
{
    if (!reverse)
    {
        scalar temp = coordinates_.b();
        coordinates_.b() = coordinates_.c();
        coordinates_.c() = coordinates_.d();
        coordinates_.d() = temp;
    }
    else
    {
        scalar temp = coordinates_.d();
        coordinates_.d() = coordinates_.c();
        coordinates_.c() = coordinates_.b();
        coordinates_.b() = temp;
    }
}


void Foam::particle::changeTet(const label tetTriI)
{
    const bool isOwner = mesh_.faceOwner()[tetFacei_] == celli_;

    const label firstTetPtI = 1;
    const label lastTetPtI = mesh_.faces()[tetFacei_].size() - 2;

    if (tetTriI == 1)
    {
        changeFace(tetTriI);
    }
    else if (tetTriI == 2)
    {
        if (isOwner)
        {
            if (tetPti_ == lastTetPtI)
            {
                changeFace(tetTriI);
            }
            else
            {
                reflect();
                tetPti_ += 1;
            }
        }
        else
        {
            if (tetPti_ == firstTetPtI)
            {
                changeFace(tetTriI);
            }
            else
            {
                reflect();
                tetPti_ -= 1;
            }
        }
    }
    else if (tetTriI == 3)
    {
        if (isOwner)
        {
            if (tetPti_ == firstTetPtI)
            {
                changeFace(tetTriI);
            }
            else
            {
                reflect();
                tetPti_ -= 1;
            }
        }
        else
        {
            if (tetPti_ == lastTetPtI)
            {
                changeFace(tetTriI);
            }
            else
            {
                reflect();
                tetPti_ += 1;
            }
        }
    }
    else
    {
        FatalErrorInFunction
            << "Changing tet without changing cell should only happen when the "
            << "track is on triangle 1, 2 or 3."
            << exit(FatalError);
    }
}


void Foam::particle::changeFace(const label tetTriI)
{
328
329
    // Get the old topology
    const triFace triOldIs(currentTetIndices().faceTriIs(mesh_));
330
331
332
333
334

    // Get the shared edge and the pre-rotation
    edge sharedEdge;
    if (tetTriI == 1)
    {
335
        sharedEdge = edge(triOldIs[1], triOldIs[2]);
336
337
338
    }
    else if (tetTriI == 2)
    {
339
        sharedEdge = edge(triOldIs[2], triOldIs[0]);
340
341
342
    }
    else if (tetTriI == 3)
    {
343
        sharedEdge = edge(triOldIs[0], triOldIs[1]);
344
345
346
347
348
349
350
    }
    else
    {
        FatalErrorInFunction
            << "Changing face without changing cell should only happen when the"
            << " track is on triangle 1, 2 or 3."
            << exit(FatalError);
351
352

        sharedEdge = edge(-1, -1);
353
354
355
356
357
358
359
360
361
    }

    // Find the face in the same cell that shares the edge, and the
    // corresponding tetrahedra point
    tetPti_ = -1;
    forAll(mesh_.cells()[celli_], cellFaceI)
    {
        const label newFaceI = mesh_.cells()[celli_][cellFaceI];
        const class face& newFace = mesh_.faces()[newFaceI];
362
        const label newOwner = mesh_.faceOwner()[newFaceI];
363
364
365
366
367
368
369

        // Exclude the current face
        if (tetFacei_ == newFaceI)
        {
            continue;
        }

370
371
372
373
        // Loop over the edges, looking for the shared one. Note that we have to
        // match the direction of the edge as well as the end points in order to
        // avoid false positives when dealing with coincident ACMI faces.
        const label edgeComp = newOwner == celli_ ? -1 : +1;
374
        label edgeI = 0;
375
376
377
378
379
380
381
        for
        (
            ;
            edgeI < newFace.size()
         && edge::compare(sharedEdge, newFace.faceEdge(edgeI)) != edgeComp;
            ++ edgeI
        );
382
383

        // If the face does not contain the edge, then move on to the next face
384
        if (edgeI >= newFace.size())
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
        {
            continue;
        }

        // Make the edge index relative to the base point
        const label newBaseI = max(0, mesh_.tetBasePtIs()[newFaceI]);
        edgeI = (edgeI - newBaseI + newFace.size()) % newFace.size();

        // If the edge is next the base point (i.e., the index is 0 or n - 1),
        // then we swap it for the adjacent edge. This new edge is opposite the
        // base point, and defines the tet with the original edge in it.
        edgeI = min(max(1, edgeI), newFace.size() - 2);

        // Set the new face and tet point
        tetFacei_ = newFaceI;
        tetPti_ = edgeI;

402
        // Exit the loop now that the tet point has been found
403
404
405
406
407
408
409
410
411
412
413
        break;
    }

    if (tetPti_ == -1)
    {
        FatalErrorInFunction
            << "The search for an edge-connected face and tet-point failed."
            << exit(FatalError);
    }

    // Pre-rotation puts the shared edge opposite the base of the tetrahedron
414
    if (sharedEdge.otherVertex(triOldIs[1]) == -1)
415
416
417
    {
        rotate(false);
    }
418
    else if (sharedEdge.otherVertex(triOldIs[2]) == -1)
419
420
421
422
    {
        rotate(true);
    }

423
424
    // Get the new topology
    const triFace triNewIs = currentTetIndices().faceTriIs(mesh_);
425
426
427
428
429

    // Reflect to account for the change of triangle orientation on the new face
    reflect();

    // Post rotation puts the shared edge back in the correct location
430
    if (sharedEdge.otherVertex(triNewIs[1]) == -1)
431
432
433
    {
        rotate(true);
    }
434
    else if (sharedEdge.otherVertex(triNewIs[2]) == -1)
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
    {
        rotate(false);
    }
}


void Foam::particle::changeCell()
{
    // Set the cell to be the one on the other side of the face
    const label ownerCellI = mesh_.faceOwner()[tetFacei_];
    const bool isOwner = celli_ == ownerCellI;
    celli_ = isOwner ? mesh_.faceNeighbour()[tetFacei_] : ownerCellI;

    // Reflect to account for the change of triangle orientation in the new cell
    reflect();
}


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
void Foam::particle::changeToMasterPatch()
{
    label thisPatch = patch();

    forAll(mesh_.cells()[celli_], cellFaceI)
    {
        // Skip the current face and any internal faces
        const label otherFaceI = mesh_.cells()[celli_][cellFaceI];
        if (facei_ == otherFaceI || mesh_.isInternalFace(otherFaceI))
        {
            continue;
        }

        // Compare the two faces. If they are the same, chose the one with the
        // lower patch index. In the case of an ACMI-wall pair, this will be
        // the ACMI, which is what we want.
        const class face& thisFace = mesh_.faces()[facei_];
        const class face& otherFace = mesh_.faces()[otherFaceI];
        if (face::compare(thisFace, otherFace) != 0)
        {
            const label otherPatch =
                mesh_.boundaryMesh().whichPatch(otherFaceI);
            if (thisPatch > otherPatch)
            {
                facei_ = otherFaceI;
                thisPatch = otherPatch;
            }
        }
    }

    tetFacei_ = facei_;
}


487
488
489
490
491
492
void Foam::particle::locate
(
    const vector& position,
    const vector* direction,
    const label celli,
    const bool boundaryFail,
493
    const string& boundaryMsg
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
)
{
    celli_ = celli;

    // Find the cell, if it has not been given
    if (celli_ < 0)
    {
        celli_ = mesh_.cellTree().findInside(position);
    }
    if (celli_ < 0)
    {
        FatalErrorInFunction
            << "Cell not found for particle position " << position << "."
            << exit(FatalError);
    }

    // Put the particle at the cell centre and in a random tet
    coordinates_ = barycentric(1, 0, 0, 0);
    tetFacei_ = mesh_.cells()[celli_][0];
    tetPti_ = 1;
    facei_ = -1;

    // Track to the injection point
    track(position - mesh_.cellCentres()[celli_], 0);
    if (!onFace())
    {
        return;
    }

    // We hit a boundary ...
    if (boundaryFail)
    {
526
527
528
529
        FatalErrorInFunction << boundaryMsg
            << " when tracking from centre " << mesh_.cellCentres()[celli_]
            << " of cell " << celli_ << " to position " << position
            << exit(FatalError);
530
531
532
533
534
535
536
537
538
    }
    else
    {
        // Re-do the track, but this time do the bit tangential to the
        // direction/patch first. This gets us as close as possible to the
        // original path/position.

        if (direction == nullptr)
        {
539
540
            const polyPatch& p = mesh_.boundaryMesh()[patch()];
            direction = &p.faceNormals()[p.whichFace(facei_)];
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
        }

        const vector n = *direction/mag(*direction);
        const vector s = position - mesh_.cellCentres()[celli_];
        const vector sN = (s & n)*n;
        const vector sT = s - sN;

        coordinates_ = barycentric(1, 0, 0, 0);
        tetFacei_ = mesh_.cells()[celli_][0];
        tetPti_ = 1;
        facei_ = -1;

        track(sT, 0);
        track(sN, 0);

        static label nWarnings = 0;
        static const label maxNWarnings = 100;
        if (nWarnings < maxNWarnings)
        {
            WarningInFunction << boundaryMsg << endl;
561
            ++nWarnings;
562
563
564
565
566
567
        }
        if (nWarnings == maxNWarnings)
        {
            WarningInFunction
                << "Suppressing any further warnings about particles being "
                << "located outside of the mesh." << endl;
568
            ++nWarnings;
569
570
        }
    }
andy's avatar
andy committed
571
572
573
574
575
576
577
578
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::particle::particle
(
    const polyMesh& mesh,
579
    const barycentric& coordinates,
580
    const label celli,
581
    const label tetFacei,
Henry Weller's avatar
Henry Weller committed
582
    const label tetPti
andy's avatar
andy committed
583
584
585
)
:
    mesh_(mesh),
586
    coordinates_(coordinates),
587
    celli_(celli),
588
    tetFacei_(tetFacei),
Henry Weller's avatar
Henry Weller committed
589
    tetPti_(tetPti),
590
591
    facei_(-1),
    stepFraction_(0.0),
andy's avatar
andy committed
592
593
594
595
596
597
598
599
600
    origProc_(Pstream::myProcNo()),
    origId_(getNewParticleID())
{}


Foam::particle::particle
(
    const polyMesh& mesh,
    const vector& position,
601
    const label celli
andy's avatar
andy committed
602
603
604
)
:
    mesh_(mesh),
605
    coordinates_(-VGREAT, -VGREAT, -VGREAT, -VGREAT),
606
    celli_(celli),
607
    tetFacei_(-1),
Henry Weller's avatar
Henry Weller committed
608
    tetPti_(-1),
609
610
    facei_(-1),
    stepFraction_(0.0),
andy's avatar
andy committed
611
612
613
    origProc_(Pstream::myProcNo()),
    origId_(getNewParticleID())
{
614
615
616
617
618
619
620
621
    locate
    (
        position,
        nullptr,
        celli,
        false,
        "Particle initialised with a location outside of the mesh."
    );
andy's avatar
andy committed
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
Foam::particle::particle
(
    const polyMesh& mesh,
    const vector& position,
    const label celli,
    const label tetFacei,
    const label tetPti,
    bool doLocate
)
:
    mesh_(mesh),
    coordinates_(-VGREAT, -VGREAT, -VGREAT, -VGREAT),
    celli_(celli),
    tetFacei_(tetFacei),
    tetPti_(tetPti),
    facei_(-1),
    stepFraction_(0.0),
    origProc_(Pstream::myProcNo()),
    origId_(getNewParticleID())
{
    if (doLocate)
    {
        locate
        (
            position,
            nullptr,
            celli,
            false,
            "Particle initialised with a location outside of the mesh."
        );
    }
}


andy's avatar
andy committed
659
660
661
Foam::particle::particle(const particle& p)
:
    mesh_(p.mesh_),
662
    coordinates_(p.coordinates_),
663
    celli_(p.celli_),
664
    tetFacei_(p.tetFacei_),
Henry Weller's avatar
Henry Weller committed
665
    tetPti_(p.tetPti_),
666
667
    facei_(p.facei_),
    stepFraction_(p.stepFraction_),
andy's avatar
andy committed
668
669
670
671
672
673
674
675
    origProc_(p.origProc_),
    origId_(p.origId_)
{}


Foam::particle::particle(const particle& p, const polyMesh& mesh)
:
    mesh_(mesh),
676
    coordinates_(p.coordinates_),
677
    celli_(p.celli_),
678
    tetFacei_(p.tetFacei_),
Henry Weller's avatar
Henry Weller committed
679
    tetPti_(p.tetPti_),
680
681
    facei_(p.facei_),
    stepFraction_(p.stepFraction_),
andy's avatar
andy committed
682
683
684
685
686
687
688
    origProc_(p.origProc_),
    origId_(p.origId_)
{}


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

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
Foam::scalar Foam::particle::track
(
    const vector& displacement,
    const scalar fraction
)
{
    scalar f = trackToFace(displacement, fraction);

    while (onInternalFace())
    {
        changeCell();

        f *= trackToFace(f*displacement, f*fraction);
    }

    return f;
}


Foam::scalar Foam::particle::trackToFace
(
    const vector& displacement,
    const scalar fraction
)
{
    scalar f = 1;

    label tetTriI = onFace() ? 0 : -1;

    facei_ = -1;

    while (true)
    {
722
        f *= trackToTri(f*displacement, f*fraction, tetTriI);
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

        if (tetTriI == -1)
        {
            // The track has completed within the current tet
            return 0;
        }
        else if (tetTriI == 0)
        {
            // The track has hit a face, so set the current face and return
            facei_ = tetFacei_;
            return f;
        }
        else
        {
            // Move to the next tet and continue the track
            changeTet(tetTriI);
        }
    }
}


Foam::scalar Foam::particle::trackToStationaryTri
(
    const vector& displacement,
    const scalar fraction,
    label& tetTriI
)
{
    const vector x0 = position();
    const vector x1 = displacement;
    const barycentric y0 = coordinates_;

    if (debug)
    {
757
758
        Info<< "Particle " << origId() << endl << "Tracking from " << x0
            << " along " << x1 << " to " << x0 + x1 << endl;
759
760
761
762
763
764
    }

    // Get the tet geometry
    vector centre;
    scalar detA;
    barycentricTensor T;
765
    stationaryTetReverseTransform(centre, detA, T);
766
767
768
769

    if (debug)
    {
        vector o, b, v1, v2;
770
        stationaryTetGeometry(o, b, v1, v2);
771
772
773
774
775
776
        Info<< "Tet points o=" << o << ", b=" << b
            << ", v1=" << v1 << ", v2=" << v2 << endl
            << "Tet determinant = " << detA << endl
            << "Start local coordinates = " << y0 << endl;
    }

777
778
779
    // Get the factor by which the displacement is increased
    const scalar f = detA >= 0 ? 1 : negativeSpaceDisplacementFactor;

780
    // Calculate the local tracking displacement
781
    barycentric Tx1(f*x1 & T);
782
783
784
785
786
787
788
789

    if (debug)
    {
        Info<< "Local displacement = " << Tx1 << "/" << detA << endl;
    }

    // Calculate the hit fraction
    label iH = -1;
790
    scalar muH = std::isnormal(detA) && detA <= 0 ? VGREAT : 1/detA;
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
    for (label i = 0; i < 4; ++ i)
    {
        if (std::isnormal(Tx1[i]) && Tx1[i] < 0)
        {
            scalar mu = - y0[i]/Tx1[i];

            if (debug)
            {
                Info<< "Hit on tet face " << i << " at local coordinate "
                    << y0 + mu*Tx1 << ", " << mu*detA*100 << "\% of the "
                    << "way along the track" << endl;
            }

            if (0 <= mu && mu < muH)
            {
                iH = i;
                muH = mu;
            }
        }
    }

    // Set the new coordinates
    barycentric yH = y0 + muH*Tx1;

    // Clamp to zero any negative coordinates generated by round-off error
    for (label i = 0; i < 4; ++ i)
    {
        yH.replace(i, i == iH ? 0 : max(0, yH[i]));
    }

    // Re-normalise if within the tet
    if (iH == -1)
    {
        yH /= cmptSum(yH);
    }

    // Set the new position and hit index
    coordinates_ = yH;
    tetTriI = iH;

    if (debug)
    {
        if (iH != -1)
        {
            Info<< "Track hit tet face " << iH << " first" << endl;
        }
        else
        {
            Info<< "Track hit no tet faces" << endl;
        }
        Info<< "End local coordinates = " << yH << endl
            << "End global coordinates = " << position() << endl
            << "Tracking displacement = " << position() - x0 << endl
            << muH*detA*100 << "\% of the step from " << stepFraction_ << " to "
            << stepFraction_ + fraction << " completed" << endl << endl;
    }

    // Set the proportion of the track that has been completed
    stepFraction_ += fraction*muH*detA;

851
    return iH != -1 ? 1 - muH*detA : 0;
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
}


Foam::scalar Foam::particle::trackToMovingTri
(
    const vector& displacement,
    const scalar fraction,
    label& tetTriI
)
{
    const vector x0 = position();
    const vector x1 = displacement;
    const barycentric y0 = coordinates_;

    if (debug)
    {
868
869
        Info<< "Particle " << origId() << endl << "Tracking from " << x0
            << " along " << x1 << " to " << x0 + x1 << endl;
870
871
872
873
874
875
876
877
    }

    // Get the tet geometry
    Pair<vector> centre;
    FixedList<scalar, 4> detA;
    FixedList<barycentricTensor, 3> T;
    movingTetReverseTransform(fraction, centre, detA, T);

878
879
880
881
    // Get the factor by which the displacement is increased
    const scalar f = detA[0] >= 0 ? 1 : negativeSpaceDisplacementFactor;

    // Get the relative global position
882
    const vector x0Rel = x0 - centre[0];
883
    const vector x1Rel = f*x1 - centre[1];
884
885
886
887

    // Form the determinant and hit equations
    cubicEqn detAEqn(sqr(detA[0])*detA[3], detA[0]*detA[2], detA[1], 1);
    const barycentric yC(1, 0, 0, 0);
888
889
    const barycentric hitEqnA =
        ((x1Rel & T[2])                  + detA[3]*yC)*sqr(detA[0]);
890
891
    const barycentric hitEqnB =
        ((x1Rel & T[1]) + (x0Rel & T[2]) + detA[2]*yC)*detA[0];
892
893
    const barycentric hitEqnC =
        ((x1Rel & T[0]) + (x0Rel & T[1]) + detA[1]*yC);
894
895
896
897
898
899
900
901
902
    const barycentric hitEqnD = y0;
    FixedList<cubicEqn, 4> hitEqn;
    forAll(hitEqn, i)
    {
        hitEqn[i] = cubicEqn(hitEqnA[i], hitEqnB[i], hitEqnC[i], hitEqnD[i]);
    }

    // Calculate the hit fraction
    label iH = -1;
903
    scalar muH = std::isnormal(detA[0]) && detA[0] <= 0 ? VGREAT : 1/detA[0];
904
    for (label i = 0; i < 4; ++i)
905
906
907
    {
        const Roots<3> mu = hitEqn[i].roots();

908
        for (label j = 0; j < 3; ++j)
909
910
911
912
913
914
915
916
917
918
919
920
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
973
974
975
976
977
978
        {
            if (mu.type(j) == roots::real && hitEqn[i].derivative(mu[j]) < 0)
            {
                if (0 <= mu[j] && mu[j] < muH)
                {
                    iH = i;
                    muH = mu[j];
                }
            }
        }
    }

    // Set the new coordinates
    barycentric yH
    (
        hitEqn[0].value(muH),
        hitEqn[1].value(muH),
        hitEqn[2].value(muH),
        hitEqn[3].value(muH)
    );
    // !!! <-- This fails if the tet collapses onto the particle, as detA tends
    // to zero at the hit. In this instance, we can differentiate the hit and
    // detA polynomials to find a limiting location, but this will not be on a
    // triangle. We will then need to do a second track through the degenerate
    // tet to find the final hit position. This second track is over zero
    // distance and therefore can be of the static mesh type. This has not yet
    // been implemented.
    const scalar detAH = detAEqn.value(muH);
    if (!std::isnormal(detAH))
    {
        FatalErrorInFunction
            << "A moving tet collapsed onto a particle. This is not supported. "
            << "The mesh is too poor, or the motion too severe, for particle "
            << "tracking to function." << exit(FatalError);
    }
    yH /= detAH;

    // Clamp to zero any negative coordinates generated by round-off error
    for (label i = 0; i < 4; ++ i)
    {
        yH.replace(i, i == iH ? 0 : max(0, yH[i]));
    }

    // Re-normalise if within the tet
    if (iH == -1)
    {
        yH /= cmptSum(yH);
    }

    // Set the new position and hit index
    coordinates_ = yH;
    tetTriI = iH;

    // Set the proportion of the track that has been completed
    stepFraction_ += fraction*muH*detA[0];

    if (debug)
    {
        if (iH != -1)
        {
            Info<< "Track hit tet face " << iH << " first" << endl;
        }
        else
        {
            Info<< "Track hit no tet faces" << endl;
        }
        Info<< "End local coordinates = " << yH << endl
            << "End global coordinates = " << position() << endl;
    }

979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
    return iH != -1 ? 1 - muH*detA[0] : 0;
}


Foam::scalar Foam::particle::trackToTri
(
    const vector& displacement,
    const scalar fraction,
    label& tetTriI
)
{
    if (mesh_.moving())
    {
        return trackToMovingTri(displacement, fraction, tetTriI);
    }
    else
    {
        return trackToStationaryTri(displacement, fraction, tetTriI);
    }
998
999
1000
}


1001
Foam::vector Foam::particle::deviationFromMeshCentre() const
1002
{
1003
    if (cmptMin(mesh_.geometricD()) == -1)
1004
    {
1005
1006
1007
        vector pos = position(), posC = pos;
        meshTools::constrainToMeshCentre(mesh_, posC);
        return pos - posC;
1008
    }
1009
    else
1010
    {
1011
        return vector::zero;
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
void Foam::particle::patchData(vector& n, vector& U) const
{
    if (!onBoundaryFace())
    {
        FatalErrorInFunction
            << "Patch data was requested for a particle that isn't on a patch"
            << exit(FatalError);
    }

    if (mesh_.moving())
    {
        Pair<vector> centre, base, vertex1, vertex2;
        movingTetGeometry(1, centre, base, vertex1, vertex2);

        n = triPointRef(base[0], vertex1[0], vertex2[0]).normal();
        n /= mag(n);

        // Interpolate the motion of the three face vertices to the current
        // coordinates
        U =
            coordinates_.b()*base[1]
          + coordinates_.c()*vertex1[1]
          + coordinates_.d()*vertex2[1];

        // The movingTetGeometry method gives the motion as a displacement
        // across the time-step, so we divide by the time-step to get velocity
        U /= mesh_.time().deltaTValue();
    }
    else
    {
        vector centre, base, vertex1, vertex2;
        stationaryTetGeometry(centre, base, vertex1, vertex2);

        n = triPointRef(base, vertex1, vertex2).normal();
        n /= mag(n);

        U = Zero;
    }
}


andy's avatar
andy committed
1057
1058
1059
1060
1061
1062
1063
1064
void Foam::particle::transformProperties(const tensor&)
{}


void Foam::particle::transformProperties(const vector&)
{}


1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
void Foam::particle::prepareForParallelTransfer()
{
    // Convert the face index to be local to the processor patch
    facei_ = mesh_.boundaryMesh()[patch()].whichFace(facei_);
}


void Foam::particle::correctAfterParallelTransfer
(
    const label patchi,
    trackingData& td
)
{
    const coupledPolyPatch& ppp =
        refCast<const coupledPolyPatch>(mesh_.boundaryMesh()[patchi]);

    if (!ppp.parallel())
    {
        const tensor& T =
        (
            ppp.forwardT().size() == 1
          ? ppp.forwardT()[0]
          : ppp.forwardT()[facei_]
        );
        transformProperties(T);
    }
    else if (ppp.separated())
    {
        const vector& s =
        (
            (ppp.separation().size() == 1)
          ? ppp.separation()[0]
          : ppp.separation()[facei_]
        );
        transformProperties(-s);
    }

    // Set the topology
    celli_ = ppp.faceCells()[facei_];
    facei_ += ppp.start();
    tetFacei_ = facei_;
    // Faces either side of a coupled patch are numbered in opposite directions
    // as their normals both point away from their connected cells. The tet
    // point therefore counts in the opposite direction from the base point.
    tetPti_ = mesh_.faces()[tetFacei_].size() - 1 - tetPti_;

    // Reflect to account for the change of triangle orientation in the new cell
    reflect();

    // Note that the position does not need transforming explicitly. The face-
    // triangle on the receive patch is the transformation of the one on the
    // send patch, so whilst the barycentric coordinates remain the same, the
    // change of triangle implicitly transforms the position.
}


1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
void Foam::particle::prepareForInteractionListReferral
(
    const vectorTensorTransform& transform
)
{
    // Get the transformed position
    const vector pos = transform.invTransformPosition(position());

    // Break the topology
    celli_ = -1;
    tetFacei_ = -1;
    tetPti_ = -1;
    facei_ = -1;

    // Store the position in the barycentric data
    coordinates_ = barycentric(1 - cmptSum(pos), pos.x(), pos.y(), pos.z());

    // Transform the properties
    transformProperties(- transform.t());
    if (transform.hasR())
    {
        transformProperties(transform.R().T());
    }
}


void Foam::particle::correctAfterInteractionListReferral(const label celli)
{
    // Get the position from the barycentric data
    const vector pos(coordinates_.b(), coordinates_.c(), coordinates_.d());

    // Create some arbitrary topology for the supplied cell
    celli_ = celli;
    tetFacei_ = mesh_.cells()[celli_][0];
    tetPti_ = 1;
    facei_ = -1;

    // Get the reverse transform and directly set the coordinates from the
    // position. This isn't likely to be correct; the particle is probably not
    // in this tet. It will, however, generate the correct vector when the
    // position method is called. A referred particle should never be tracked,
    // so this approximate topology is good enough. By using the nearby cell we
    // minimise the error associated with the incorrect topology.
1164
    coordinates_ = barycentric(1, 0, 0, 0);
1165
1166
    if (mesh_.moving())
    {
1167
1168
1169
1170
1171
        Pair<vector> centre;
        FixedList<scalar, 4> detA;
        FixedList<barycentricTensor, 3> T;
        movingTetReverseTransform(0, centre, detA, T);
        coordinates_ += (pos - centre[0]) & T[0]/detA[0];
1172
1173
1174
1175
1176
1177
    }
    else
    {
        vector centre;
        scalar detA;
        barycentricTensor T;
1178
1179
        stationaryTetReverseTransform(centre, detA, T);
        coordinates_ += (pos - centre) & T/detA;
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
    }
}


Foam::label Foam::particle::procTetPt
(
    const polyMesh& procMesh,
    const label procCell,
    const label procTetFace
) const
{
    // The tet point on the procMesh differs from the current tet point if the
    // mesh and procMesh faces are of differing orientation. The change is the
    // same as in particle::correctAfterParallelTransfer.

    if
    (
        (mesh_.faceOwner()[tetFacei_] == celli_)
     == (procMesh.faceOwner()[procTetFace] == procCell)
    )
    {
        return tetPti_;
    }
    else
    {
        return procMesh.faces()[procTetFace].size() - 1 - tetPti_;
    }
}


void Foam::particle::autoMap
(
    const vector& position,
    const mapPolyMesh& mapper
)
{
    locate
    (
        position,
        nullptr,
        mapper.reverseCellMap()[celli_],
        true,
        "Particle mapped to a location outside of the mesh."
    );
}


1227
void Foam::particle::relocate(const point& position)
1228
1229
1230
{
    locate
    (
1231
        position,
1232
1233
1234
1235
1236
1237
1238
1239
        nullptr,
        celli_,
        true,
        "Particle mapped to a location outside of the mesh."
    );
}


andy's avatar
andy committed
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
// * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * * //

bool Foam::operator==(const particle& pA, const particle& pB)
{
    return (pA.origProc() == pB.origProc() && pA.origId() == pB.origId());
}


bool Foam::operator!=(const particle& pA, const particle& pB)
{
    return !(pA == pB);
}


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