sampledSet.C 11.8 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2015 OpenFOAM Foundation
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 "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 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59

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

namespace Foam
{
    const scalar sampledSet::tol = 1e-6;

    defineTypeNameAndDebug(sampledSet, 0);
    defineRunTimeSelectionTable(sampledSet, word);
}


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

Foam::label Foam::sampledSet::getBoundaryCell(const label faceI) const
{
    return mesh().faceOwner()[faceI];
}


Foam::label Foam::sampledSet::getCell
(
    const label faceI,
    const point& sample
) const
{
    if (faceI == -1)
    {
60 61
        FatalErrorInFunction
            << "Illegal face label " << faceI
62 63 64 65 66 67 68
            << abort(FatalError);
    }

    if (faceI >= mesh().nInternalFaces())
    {
        label cellI = getBoundaryCell(faceI);

69
        if (!mesh().pointInCell(sample, cellI, searchEngine_.decompMode()))
70
        {
71 72
            FatalErrorInFunction
                << "Found cell " << cellI << " using face " << faceI
73 74 75 76 77 78 79 80 81 82 83
                << ". But cell does not contain point " << sample
                << abort(FatalError);
        }
        return cellI;
    }
    else
    {
        // Try owner and neighbour to see which one contains sample

        label cellI = mesh().faceOwner()[faceI];

84
        if (mesh().pointInCell(sample, cellI, searchEngine_.decompMode()))
85 86 87 88 89 90 91
        {
            return cellI;
        }
        else
        {
            cellI = mesh().faceNeighbour()[faceI];

92
            if (mesh().pointInCell(sample, cellI, searchEngine_.decompMode()))
93 94 95 96 97
            {
                return cellI;
            }
            else
            {
98 99
                FatalErrorInFunction
                    << "None of the neighbours of face "
100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180
                    << faceI << " contains point " << sample
                    << abort(FatalError);

                return -1;
            }
        }
    }
}


Foam::scalar Foam::sampledSet::calcSign
(
    const label faceI,
    const point& sample
) const
{
    vector vec = sample - mesh().faceCentres()[faceI];

    scalar magVec = mag(vec);

    if (magVec < VSMALL)
    {
        // sample on face centre. Regard as inside
        return -1;
    }

    vec /= magVec;

    vector n = mesh().faceAreas()[faceI];

    n /= mag(n) + VSMALL;

    return n & vec;
}


// Return face (or -1) of face which is within smallDist of sample
Foam::label Foam::sampledSet::findNearFace
(
    const label cellI,
    const point& sample,
    const scalar smallDist
) const
{
    const cell& myFaces = mesh().cells()[cellI];

    forAll(myFaces, myFaceI)
    {
        const face& f = mesh().faces()[myFaces[myFaceI]];

        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)
        {
            return myFaces[myFaceI];
        }
    }
    return -1;
}


// 'Pushes' point facePt (which is almost on face) in direction of cell centre
// so it is clearly inside.
Foam::point Foam::sampledSet::pushIn
(
    const point& facePt,
    const label faceI
) const
{
    label cellI = mesh().faceOwner()[faceI];
181
    const point& cC = mesh().cellCentres()[cellI];
182

183
    point newPosition = facePt;
184

185 186 187 188
    // Taken from particle::initCellFacePt()
    label tetFaceI;
    label tetPtI;
    mesh().findTetFacePt(cellI, facePt, tetFaceI, tetPtI);
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
    if (tetFaceI == -1 || tetPtI == -1)
    {
        newPosition = facePt;

        label trap(1.0/particle::trackingCorrectionTol + 1);

        label iterNo = 0;

        do
        {
            newPosition += particle::trackingCorrectionTol*(cC - facePt);

            mesh().findTetFacePt
            (
                cellI,
                newPosition,
                tetFaceI,
                tetPtI
            );

            iterNo++;

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

    if (tetFaceI == -1)
216
    {
217 218
        FatalErrorInFunction
            << "After pushing " << facePt << " to " << newPosition
219 220 221 222
            << " it is still outside face " << faceI
            << " at " << mesh().faceCentres()[faceI]
            << " of cell " << cellI
            << " at " << cC << endl
223 224 225
            << "Please change your starting point"
            << abort(FatalError);
    }
226 227

    //Info<< "pushIn : moved " << facePt << " to " << newPosition
228 229
    //    << endl;

230
    return newPosition;
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
}


// Calculates start of tracking given samplePt and first boundary intersection
// (bPoint, bFaceI). bFaceI == -1 if no boundary intersection.
// Returns true if trackPt is sampling point
bool Foam::sampledSet::getTrackingPoint
(
    const vector& offset,
    const point& samplePt,
    const point& bPoint,
    const label bFaceI,

    point& trackPt,
    label& trackCellI,
    label& trackFaceI
) const
{
    const scalar smallDist = mag(tol*offset);

    bool isGoodSample = false;

    if (bFaceI == -1)
    {
        // No boundary intersection. Try and find cell samplePt is in
256
        trackCellI = mesh().findCell(samplePt, searchEngine_.decompMode());
257 258 259 260

        if
        (
            (trackCellI == -1)
261 262 263 264 265 266
        || !mesh().pointInCell
            (
                samplePt,
                trackCellI,
                searchEngine_.decompMode()
            )
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
        )
        {
            // Line samplePt - end_ does not intersect domain at all.
            // (or is along edge)
            //Info<< "getTrackingPoint : samplePt outside domain : "
            //    << "  samplePt:" << samplePt
            //    << endl;

            trackCellI = -1;
            trackFaceI = -1;

            isGoodSample = false;
        }
        else
        {
            // start is inside. Use it as tracking point
            //Info<< "getTrackingPoint : samplePt inside :"
            //    << "  samplePt:" << samplePt
            //    << "  trackCellI:" << trackCellI
            //    << endl;

            trackPt = samplePt;
            trackFaceI = -1;

            isGoodSample = true;
        }
    }
    else if (mag(samplePt - bPoint) < smallDist)
    {
        //Info<< "getTrackingPoint : samplePt:" << samplePt
        //    << " close to bPoint:"
        //    << bPoint << endl;

        // samplePt close to bPoint. Snap to it
        trackPt = pushIn(bPoint, bFaceI);
        trackFaceI = bFaceI;
        trackCellI = getBoundaryCell(trackFaceI);

        isGoodSample = true;
    }
    else
    {
        scalar sign = calcSign(bFaceI, samplePt);

        if (sign < 0)
        {
            // samplePt inside or marginally outside.
            trackPt = samplePt;
            trackFaceI = -1;
316
            trackCellI = mesh().findCell(trackPt, searchEngine_.decompMode());
317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372

            isGoodSample = true;
        }
        else
        {
            // samplePt outside. use bPoint
            trackPt = pushIn(bPoint, bFaceI);
            trackFaceI = bFaceI;
            trackCellI = getBoundaryCell(trackFaceI);

            isGoodSample = false;
        }
    }

    if (debug)
    {
        Info<< "sampledSet::getTrackingPoint :"
            << " offset:" << offset
            << " samplePt:" << samplePt
            << " bPoint:" << bPoint
            << " bFaceI:" << bFaceI
            << endl << "   Calculated first tracking point :"
            << " trackPt:" << trackPt
            << " trackCellI:" << trackCellI
            << " trackFaceI:" << trackFaceI
            << " 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
)
{
    setSize(samplingPts.size());
    cells_.setSize(samplingCells.size());
    faces_.setSize(samplingFaces.size());
    segments_.setSize(samplingSegments.size());
    curveDist_.setSize(samplingCurveDist.size());

    if
    (
        (cells_.size() != size())
     || (faces_.size() != size())
     || (segments_.size() != size())
     || (curveDist_.size() != size())
    )
    {
373
        FatalErrorInFunction
374 375 376 377 378 379 380 381 382 383 384 385 386
            << "sizes not equal : "
            << "  points:" << size()
            << "  cells:" << cells_.size()
            << "  faces:" << faces_.size()
            << "  segments:" << segments_.size()
            << "  curveDist:" << curveDist_.size()
            << abort(FatalError);
    }

    forAll(samplingPts, sampleI)
    {
        operator[](sampleI) = samplingPts[sampleI];
    }
387 388
    curveDist_ = samplingCurveDist;

389 390 391 392 393 394 395 396 397 398 399 400
    cells_ = samplingCells;
    faces_ = samplingFaces;
    segments_ = samplingSegments;
}


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

Foam::sampledSet::sampledSet
(
    const word& name,
    const polyMesh& mesh,
401
    const meshSearch& searchEngine,
402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417
    const word& axis
)
:
    coordSet(name, axis),
    mesh_(mesh),
    searchEngine_(searchEngine),
    segments_(0),
    cells_(0),
    faces_(0)
{}


Foam::sampledSet::sampledSet
(
    const word& name,
    const polyMesh& mesh,
418
    const meshSearch& searchEngine,
419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438
    const dictionary& dict
)
:
    coordSet(name, dict.lookup("axis")),
    mesh_(mesh),
    searchEngine_(searchEngine),
    segments_(0),
    cells_(0),
    faces_(0)
{}


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

Foam::sampledSet::~sampledSet()
{}


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

439
Foam::autoPtr<Foam::sampledSet> Foam::sampledSet::New
440 441 442
(
    const word& name,
    const polyMesh& mesh,
443
    const meshSearch& searchEngine,
444 445 446
    const dictionary& dict
)
{
447
    const word sampleType(dict.lookup("type"));
448 449 450 451 452 453

    wordConstructorTable::iterator cstrIter =
        wordConstructorTablePtr_->find(sampleType);

    if (cstrIter == wordConstructorTablePtr_->end())
    {
454 455
        FatalErrorInFunction
            << "Unknown sample type "
456
            << sampleType << nl << nl
457
            << "Valid sample types : " << endl
458
            << wordConstructorTablePtr_->sortedToc()
459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492
            << exit(FatalError);
    }

    return autoPtr<sampledSet>
    (
        cstrIter()
        (
            name,
            mesh,
            searchEngine,
            dict
        )
    );
}


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

    os  << endl << "\t(cellI)\t(faceI)" << endl;

    forAll(*this, sampleI)
    {
        os  << '\t' << cells_[sampleI]
            << '\t' << faces_[sampleI]
            << endl;
    }

    return os;
}


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