sampledSet.C 14.2 KB
Newer Older
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.
7
8
9
10
-------------------------------------------------------------------------------
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 "sampledSet.H"
#include "polyMesh.H"
#include "primitiveMesh.H"
#include "meshSearch.H"
#include "writer.H"
31
#include "particle.H"
32
33
34
35
36
37
38
39
40
41

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

namespace Foam
{
    defineTypeNameAndDebug(sampledSet, 0);
    defineRunTimeSelectionTable(sampledSet, word);
}


Henry Weller's avatar
Henry Weller committed
42
// * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
43

44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
void Foam::sampledSet::checkDimensions() const
{
    if
    (
        (cells_.size() != size())
     || (faces_.size() != size())
     || (segments_.size() != size())
     || (curveDist_.size() != size())
    )
    {
        FatalErrorInFunction
            << "sizes not equal : "
            << "  points:" << size()
            << "  cells:" << cells_.size()
            << "  faces:" << faces_.size()
            << "  segments:" << segments_.size()
            << "  curveDist:" << curveDist_.size()
            << abort(FatalError);
    }
}


66
Foam::label Foam::sampledSet::getBoundaryCell(const label facei) const
67
{
68
    return mesh().faceOwner()[facei];
69
70
71
}


72
Foam::label Foam::sampledSet::getNeighbourCell(const label facei) const
73
{
74
    if (facei >= mesh().nInternalFaces())
75
    {
76
        return mesh().faceOwner()[facei];
77
78
79
    }
    else
    {
80
        return mesh().faceNeighbour()[facei];
81
82
83
84
    }
}


85
Foam::label Foam::sampledSet::pointInCell
86
(
87
88
    const point& p,
    const label samplei
89
90
) const
{
91
92
    // Collect the face owner and neighbour cells of the sample into an array
    // for convenience
93
    const label cells[4] =
94
95
96
97
98
99
100
101
102
    {
        mesh().faceOwner()[faces_[samplei]],
        getNeighbourCell(faces_[samplei]),
        mesh().faceOwner()[faces_[samplei+1]],
        getNeighbourCell(faces_[samplei+1])
    };

    // Find the sampled cell by checking the owners and neighbours of the
    // sampled faces
103
    label cellm =
104
105
        (cells[0] == cells[2] || cells[0] == cells[3]) ? cells[0]
      : (cells[1] == cells[2] || cells[1] == cells[3]) ? cells[1]
106
      : -1;
107

108
109
    if (cellm != -1)
    {
110
111
        // If found the sampled cell check the point is in the cell
        // otherwise ignore
112
        if (!mesh().pointInCell(p, cellm, searchEngine_.decompMode()))
113
        {
114
            cellm = -1;
115

116
            if (debug)
117
            {
118
119
120
                WarningInFunction
                    << "Could not find mid-point " << p
                    << " cell " << cellm << endl;
121
122
123
            }
        }
    }
124
    else
125
    {
126
127
        // If the sample does not pass through a single cell check if the point
        // is in any of the owners or neighbours otherwise ignore
128
        for (label i=0; i<4; ++i)
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
        {
            if (mesh().pointInCell(p, cells[i], searchEngine_.decompMode()))
            {
                return cells[i];
            }
        }

        if (debug)
        {
            WarningInFunction
                << "Could not find cell for mid-point" << nl
                << "  samplei: " << samplei
                << "  pts[samplei]: " << operator[](samplei)
                << "  face[samplei]: " << faces_[samplei]
                << "  pts[samplei+1]: " << operator[](samplei+1)
                << "  face[samplei+1]: " << faces_[samplei+1]
                << "  cellio: " << cells[0]
                << "  cellin: " << cells[1]
                << "  celljo: " << cells[2]
                << "  celljn: " << cells[3]
                << endl;
        }
151
152
153
    }

    return cellm;
154
155
156
157
158
}


Foam::scalar Foam::sampledSet::calcSign
(
159
    const label facei,
160
161
162
    const point& sample
) const
{
163
    vector vec = sample - mesh().faceCentres()[facei];
164
165
166
167
168

    scalar magVec = mag(vec);

    if (magVec < VSMALL)
    {
169
        // Sample on face centre. Regard as inside
170
171
172
173
174
        return -1;
    }

    vec /= magVec;

175
    vector n = mesh().faceAreas()[facei];
176
177
178
179
180
181
182
183
184

    n /= mag(n) + VSMALL;

    return n & vec;
}


Foam::label Foam::sampledSet::findNearFace
(
185
    const label celli,
186
187
188
189
    const point& sample,
    const scalar smallDist
) const
{
190
    const cell& myFaces = mesh().cells()[celli];
191

192
    forAll(myFaces, myFacei)
193
    {
194
        const face& f = mesh().faces()[myFaces[myFacei]];
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210

        pointHit inter = f.nearestPoint(sample, mesh().points());

        scalar dist;

        if (inter.hit())
        {
            dist = mag(inter.hitPoint() - sample);
        }
        else
        {
            dist = mag(inter.missPoint() - sample);
        }

        if (dist < smallDist)
        {
211
            return myFaces[myFacei];
212
213
214
215
216
217
218
219
220
        }
    }
    return -1;
}


Foam::point Foam::sampledSet::pushIn
(
    const point& facePt,
221
    const label facei
222
223
) const
{
224
225
    label celli = mesh().faceOwner()[facei];
    const point& cC = mesh().cellCentres()[celli];
226

227
    point newPosition = facePt;
228

229
    // Taken from particle::initCellFacePt()
230
    label tetFacei;
231
    label tetPtI;
232
    mesh().findTetFacePt(celli, facePt, tetFacei, tetPtI);
233

234
235
236
237
238
239
240
241
    // This is the tolerance that was defined as a static constant of the
    // particle class. It is no longer used by particle, following the switch to
    // barycentric tracking. The only place in which the tolerance is now used
    // is here. I'm not sure what the purpose of this code is, but it probably
    // wants removing. It is doing tet-searches for a particle position. This
    // should almost certainly be left to the particle class.
    const scalar trackingCorrectionTol = 1e-5;

242
    if (tetFacei == -1 || tetPtI == -1)
243
244
245
    {
        newPosition = facePt;

246
        label trap(1.0/trackingCorrectionTol + 1);
247
248
249
250
251

        label iterNo = 0;

        do
        {
252
            newPosition += trackingCorrectionTol*(cC - facePt);
253
254
255

            mesh().findTetFacePt
            (
256
                celli,
257
                newPosition,
258
                tetFacei,
259
260
261
                tetPtI
            );

262
            ++iterNo;
263

264
        } while (tetFacei < 0  && iterNo <= trap);
265
266
    }

267
    if (tetFacei == -1)
268
    {
269
270
        FatalErrorInFunction
            << "After pushing " << facePt << " to " << newPosition
271
272
273
            << " it is still outside face " << facei
            << " at " << mesh().faceCentres()[facei]
            << " of cell " << celli
274
            << " at " << cC << endl
275
276
277
            << "Please change your starting point"
            << abort(FatalError);
    }
278
279

    return newPosition;
280
281
282
283
284
285
286
}


bool Foam::sampledSet::getTrackingPoint
(
    const point& samplePt,
    const point& bPoint,
287
    const label bFacei,
288
    const scalar smallDist,
289
290

    point& trackPt,
291
292
    label& trackCelli,
    label& trackFacei
293
294
295
296
) const
{
    bool isGoodSample = false;

297
    if (bFacei == -1)
298
299
    {
        // No boundary intersection. Try and find cell samplePt is in
300
        trackCelli = mesh().findCell(samplePt, searchEngine_.decompMode());
301
302
303

        if
        (
304
            (trackCelli == -1)
305
306
307
        || !mesh().pointInCell
            (
                samplePt,
308
                trackCelli,
309
310
                searchEngine_.decompMode()
            )
311
312
313
314
315
        )
        {
            // Line samplePt - end_ does not intersect domain at all.
            // (or is along edge)

316
317
            trackCelli = -1;
            trackFacei = -1;
318
319
320
321
322

            isGoodSample = false;
        }
        else
        {
323
            // Start is inside. Use it as tracking point
324
325

            trackPt = samplePt;
326
            trackFacei = -1;
327
328
329
330
331
332
333

            isGoodSample = true;
        }
    }
    else if (mag(samplePt - bPoint) < smallDist)
    {
        // samplePt close to bPoint. Snap to it
334
335
336
        trackPt = pushIn(bPoint, bFacei);
        trackFacei = bFacei;
        trackCelli = getBoundaryCell(trackFacei);
337
338
339
340
341

        isGoodSample = true;
    }
    else
    {
342
        scalar sign = calcSign(bFacei, samplePt);
343
344
345
346
347

        if (sign < 0)
        {
            // samplePt inside or marginally outside.
            trackPt = samplePt;
348
349
            trackFacei = -1;
            trackCelli = mesh().findCell(trackPt, searchEngine_.decompMode());
350
351
352
353
354
355

            isGoodSample = true;
        }
        else
        {
            // samplePt outside. use bPoint
356
357
358
            trackPt = pushIn(bPoint, bFacei);
            trackFacei = bFacei;
            trackCelli = getBoundaryCell(trackFacei);
359
360
361
362
363
364
365

            isGoodSample = false;
        }
    }

    if (debug)
    {
366
        InfoInFunction
367
368
            << " samplePt:" << samplePt
            << " bPoint:" << bPoint
369
            << " bFacei:" << bFacei
370
371
            << endl << "   Calculated first tracking point :"
            << " trackPt:" << trackPt
372
373
            << " trackCelli:" << trackCelli
            << " trackFacei:" << trackFacei
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
            << " isGoodSample:" << isGoodSample
            << endl;
    }

    return isGoodSample;
}


void Foam::sampledSet::setSamples
(
    const List<point>& samplingPts,
    const labelList& samplingCells,
    const labelList& samplingFaces,
    const labelList& samplingSegments,
    const scalarList& samplingCurveDist
)
{
391
    setPoints(samplingPts);
392
393
    curveDist_ = samplingCurveDist;

394
    segments_ = samplingSegments;
395
396
    cells_ = samplingCells;
    faces_ = samplingFaces;
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418

    checkDimensions();
}


void Foam::sampledSet::setSamples
(
    List<point>&& samplingPts,
    labelList&& samplingCells,
    labelList&& samplingFaces,
    labelList&& samplingSegments,
    scalarList&& samplingCurveDist
)
{
    setPoints(std::move(samplingPts));
    curveDist_ = std::move(samplingCurveDist);

    segments_ = std::move(samplingSegments);
    cells_ = std::move(samplingCells);
    faces_ = std::move(samplingFaces);

    checkDimensions();
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
Foam::autoPtr<Foam::coordSet> Foam::sampledSet::gather
(
    labelList& indexSet
) const
{
    // Combine sampleSet from processors. Sort by curveDist. Return
    // ordering in indexSet.
    // Note: only master results are valid

    // Collect data from all processors
    List<List<point>> gatheredPts(Pstream::nProcs());
    gatheredPts[Pstream::myProcNo()] = *this;
    Pstream::gatherList(gatheredPts);

    List<labelList> gatheredSegments(Pstream::nProcs());
    gatheredSegments[Pstream::myProcNo()] = segments();
    Pstream::gatherList(gatheredSegments);

    List<scalarList> gatheredDist(Pstream::nProcs());
    gatheredDist[Pstream::myProcNo()] = curveDist();
    Pstream::gatherList(gatheredDist);


    // Combine processor lists into one big list.
    List<point> allPts
    (
        ListListOps::combine<List<point>>
        (
            gatheredPts, accessOp<List<point>>()
        )
    );
    labelList allSegments
    (
        ListListOps::combine<labelList>
        (
            gatheredSegments, accessOp<labelList>()
        )
    );
    scalarList allCurveDist
    (
        ListListOps::combine<scalarList>
        (
            gatheredDist, accessOp<scalarList>()
        )
    );


469
    if (Pstream::master() && allCurveDist.empty())
470
471
472
473
474
475
476
    {
        WarningInFunction
            << "Sample set " << name()
            << " has zero points." << endl;
    }

    // Sort curveDist and use to fill masterSamplePts
477
478
    Foam::sortedOrder(allCurveDist, indexSet);      // uses stable sort
    scalarList sortedDist(allCurveDist, indexSet);  // with indices for mapping
479

480
    return autoPtr<coordSet>::New
481
    (
482
483
484
485
        name(),
        axis(),
        List<point>(UIndirectList<point>(allPts, indexSet)),
        sortedDist
486
487
488
489
    );
}


490
491
492
493
494
495
// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::sampledSet::sampledSet
(
    const word& name,
    const polyMesh& mesh,
496
    const meshSearch& searchEngine,
497
    const coordSet::coordFormat axisType
498
499
)
:
500
    coordSet(name, axisType),
501
502
    mesh_(mesh),
    searchEngine_(searchEngine),
503
504
505
    segments_(),
    cells_(),
    faces_()
506
507
508
509
510
511
512
{}


Foam::sampledSet::sampledSet
(
    const word& name,
    const polyMesh& mesh,
513
    const meshSearch& searchEngine,
514
    const word& axis
515
516
)
:
517
    coordSet(name, axis),
518
519
    mesh_(mesh),
    searchEngine_(searchEngine),
520
521
522
    segments_(),
    cells_(),
    faces_()
523
524
525
{}


526
527
528
529
530
531
532
533
534
535
536
537
538
539
Foam::sampledSet::sampledSet
(
    const word& name,
    const polyMesh& mesh,
    const meshSearch& searchEngine,
    const dictionary& dict
)
:
    coordSet(name, dict.get<word>("axis")),
    mesh_(mesh),
    searchEngine_(searchEngine),
    segments_(),
    cells_(),
    faces_()
540
541
542
543
544
{}


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

545
Foam::autoPtr<Foam::sampledSet> Foam::sampledSet::New
546
547
548
(
    const word& name,
    const polyMesh& mesh,
549
    const meshSearch& searchEngine,
550
551
552
    const dictionary& dict
)
{
553
    const word sampleType(dict.get<word>("type"));
554

555
    auto cstrIter = wordConstructorTablePtr_->cfind(sampleType);
556

557
    if (!cstrIter.found())
558
    {
559
560
        FatalErrorInFunction
            << "Unknown sample type "
561
            << sampleType << nl << nl
562
            << "Valid sample types : " << nl
563
            << wordConstructorTablePtr_->sortedToc()
564
565
566
567
568
569
570
571
572
573
            << exit(FatalError);
    }

    return autoPtr<sampledSet>
    (
        cstrIter()
        (
            name,
            mesh,
            searchEngine,
574
            dict.optionalSubDict(sampleType + "Coeffs")
575
576
577
578
579
580
581
582
583
        )
    );
}


Foam::Ostream& Foam::sampledSet::write(Ostream& os) const
{
    coordSet::write(os);

584
    os  << nl << "\t(celli)\t(facei)" << nl;
585

586
    forAll(*this, samplei)
587
    {
588
589
590
        os  << '\t' << cells_[samplei]
            << '\t' << faces_[samplei]
            << nl;
591
592
593
594
595
596
597
    }

    return os;
}


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