Skip to content
Snippets Groups Projects
MRFZone.C 14.4 KiB
Newer Older
  • Learn to ignore specific revisions
  • /*---------------------------------------------------------------------------*\
      =========                 |
      \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
       \\    /   O peration     |
    
    OpenFOAM bot's avatar
    OpenFOAM bot committed
        \\  /    A nd           | www.openfoam.com
    
         \\/     M anipulation  |
    -------------------------------------------------------------------------------
    
    OpenFOAM bot's avatar
    OpenFOAM bot committed
        Copyright (C) 2011-2017 OpenFOAM Foundation
    
        Copyright (C) 2018-2021 OpenCFD Ltd.
    
    -------------------------------------------------------------------------------
    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 "MRFZone.H"
    #include "fvMesh.H"
    #include "volFields.H"
    #include "surfaceFields.H"
    #include "fvMatrices.H"
    
    mattijs's avatar
    mattijs committed
    #include "faceSet.H"
    
    #include "geometricOneField.H"
    
    #include "syncTools.H"
    
    mattijs's avatar
    mattijs committed
    
    
    // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
    
    
    namespace Foam
    {
    
    mattijs's avatar
    mattijs committed
    // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
    
    
    void Foam::MRFZone::setMRFFaces()
    
    mattijs's avatar
    mattijs committed
        const polyBoundaryMesh& patches = mesh_.boundaryMesh();
    
    
        // Type per face:
        //  0:not in zone
        //  1:moving with frame
        //  2:other
    
        labelList faceType(mesh_.nFaces(), Zero);
    
    
        // Determine faces in cell zone
        // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
        // (without constructing cells)
    
        const labelList& own = mesh_.faceOwner();
        const labelList& nei = mesh_.faceNeighbour();
    
        // Cells in zone
        boolList zoneCell(mesh_.nCells(), false);
    
        if (cellZoneID_ != -1)
    
    mattijs's avatar
    mattijs committed
        {
    
            const labelList& cellLabels = mesh_.cellZones()[cellZoneID_];
            forAll(cellLabels, i)
    
    mattijs's avatar
    mattijs committed
            {
    
                zoneCell[cellLabels[i]] = true;
            }
        }
    
        label nZoneFaces = 0;
    
    
        for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
    
            if (zoneCell[own[facei]] || zoneCell[nei[facei]])
    
                faceType[facei] = 1;
    
    
        labelHashSet excludedPatches(excludedPatchLabels_);
    
    
        forAll(patches, patchi)
    
            const polyPatch& pp = patches[patchi];
    
    mattijs's avatar
    mattijs committed
    
    
            if (pp.coupled() || excludedPatches.found(patchi))
    
    mattijs's avatar
    mattijs committed
            {
    
                    label facei = pp.start()+i;
    
    mattijs's avatar
    mattijs committed
    
    
                    if (zoneCell[own[facei]])
    
                        faceType[facei] = 2;
    
                        nZoneFaces++;
                    }
                }
            }
            else if (!isA<emptyPolyPatch>(pp))
            {
                forAll(pp, i)
    
    mattijs's avatar
    mattijs committed
                {
    
                    label facei = pp.start()+i;
    
                    if (zoneCell[own[facei]])
    
                        faceType[facei] = 1;
    
    mattijs's avatar
    mattijs committed
                }
            }
    
        // Synchronize the faceType across processor patches
        syncTools::syncFaceList(mesh_, faceType, maxEqOp<label>());
    
    
        // Now we have for faceType:
        //  0   : face not in cellZone
        //  1   : internal face or normal patch face
        //  2   : coupled patch face or excluded patch face
    
        // Sort into lists per patch.
    
    
    mattijs's avatar
    mattijs committed
        internalFaces_.setSize(mesh_.nFaces());
        label nInternal = 0;
    
    
        for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
    
    mattijs's avatar
    mattijs committed
        {
    
            if (faceType[facei] == 1)
    
    mattijs's avatar
    mattijs committed
            {
    
                internalFaces_[nInternal++] = facei;
    
    mattijs's avatar
    mattijs committed
            }
        }
        internalFaces_.setSize(nInternal);
    
    
        labelList nIncludedFaces(patches.size(), Zero);
        labelList nExcludedFaces(patches.size(), Zero);
    
    mattijs's avatar
    mattijs committed
    
        forAll(patches, patchi)
        {
            const polyPatch& pp = patches[patchi];
    
            forAll(pp, patchFacei)
            {
    
                label facei = pp.start() + patchFacei;
    
    mattijs's avatar
    mattijs committed
    
    
                if (faceType[facei] == 1)
    
    mattijs's avatar
    mattijs committed
                {
                    nIncludedFaces[patchi]++;
                }
    
                else if (faceType[facei] == 2)
    
    mattijs's avatar
    mattijs committed
                {
                    nExcludedFaces[patchi]++;
                }
            }
        }
    
        includedFaces_.setSize(patches.size());
        excludedFaces_.setSize(patches.size());
        forAll(nIncludedFaces, patchi)
        {
            includedFaces_[patchi].setSize(nIncludedFaces[patchi]);
            excludedFaces_[patchi].setSize(nExcludedFaces[patchi]);
        }
        nIncludedFaces = 0;
        nExcludedFaces = 0;
    
        forAll(patches, patchi)
        {
            const polyPatch& pp = patches[patchi];
    
            forAll(pp, patchFacei)
            {
    
                label facei = pp.start() + patchFacei;
    
    mattijs's avatar
    mattijs committed
    
    
                if (faceType[facei] == 1)
    
    mattijs's avatar
    mattijs committed
                {
                    includedFaces_[patchi][nIncludedFaces[patchi]++] = patchFacei;
                }
    
                else if (faceType[facei] == 2)
    
    mattijs's avatar
    mattijs committed
                {
                    excludedFaces_[patchi][nExcludedFaces[patchi]++] = patchFacei;
                }
            }
        }
    
    
    
            faceSet internalFaces(mesh_, "internalFaces", internalFaces_);
            Pout<< "Writing " << internalFaces.size()
                << " internal faces in MRF zone to faceSet "
                << internalFaces.name() << endl;
            internalFaces.write();
    
            faceSet MRFFaces(mesh_, "includedFaces", 100);
            forAll(includedFaces_, patchi)
    
    mattijs's avatar
    mattijs committed
            {
    
                forAll(includedFaces_[patchi], i)
    
    mattijs's avatar
    mattijs committed
                {
    
                    label patchFacei = includedFaces_[patchi][i];
                    MRFFaces.insert(patches[patchi].start()+patchFacei);
    
    mattijs's avatar
    mattijs committed
                }
            }
    
            Pout<< "Writing " << MRFFaces.size()
                << " patch faces in MRF zone to faceSet "
                << MRFFaces.name() << endl;
            MRFFaces.write();
    
    mattijs's avatar
    mattijs committed
    
    
            faceSet excludedFaces(mesh_, "excludedFaces", 100);
            forAll(excludedFaces_, patchi)
    
    mattijs's avatar
    mattijs committed
            {
    
                forAll(excludedFaces_[patchi], i)
    
    mattijs's avatar
    mattijs committed
                {
    
                    label patchFacei = excludedFaces_[patchi][i];
                    excludedFaces.insert(patches[patchi].start()+patchFacei);
    
    mattijs's avatar
    mattijs committed
                }
            }
    
            Pout<< "Writing " << excludedFaces.size()
                << " faces in MRF zone with special handling to faceSet "
                << excludedFaces.name() << endl;
            excludedFaces.write();
        }
    
    mattijs's avatar
    mattijs committed
    }
    
    mattijs's avatar
    mattijs committed
    
    
    
    mattijs's avatar
    mattijs committed
    // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
    
    mattijs's avatar
    mattijs committed
    
    
    andy's avatar
    andy committed
    Foam::MRFZone::MRFZone
    (
        const word& name,
        const fvMesh& mesh,
    
    andy's avatar
    andy committed
        const dictionary& dict,
        const word& cellZoneName
    
    andy's avatar
    andy committed
    )
    
    mattijs's avatar
    mattijs committed
    :
        mesh_(mesh),
    
    andy's avatar
    andy committed
        name_(name),
        coeffs_(dict),
    
    andy's avatar
    andy committed
        cellZoneName_(cellZoneName),
    
        cellZoneID_(-1),
        excludedPatchNames_(wordRes()),
        origin_(Zero),
        axis_(Zero),
        omega_(nullptr)
    
    mattijs's avatar
    mattijs committed
    {
    
    }
    
    
    // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
    
    
    Foam::vector Foam::MRFZone::Omega() const
    {
        return omega_->value(mesh_.time().timeOutputValue())*axis_;
    }
    
    
    
    void Foam::MRFZone::addCoriolis
    (
        const volVectorField& U,
        volVectorField& ddtU
    ) const
    {
        if (cellZoneID_ == -1)
        {
            return;
        }
    
        const labelList& cells = mesh_.cellZones()[cellZoneID_];
    
        vectorField& ddtUc = ddtU.primitiveFieldRef();
    
        const vector Omega = this->Omega();
    
    
        forAll(cells, i)
        {
            label celli = cells[i];
    
            ddtUc[celli] += (Omega ^ Uc[celli]);
    
    andy's avatar
    andy committed
    void Foam::MRFZone::addCoriolis(fvVectorMatrix& UEqn, const bool rhs) const
    
    {
        if (cellZoneID_ == -1)
        {
            return;
        }
    
        const labelList& cells = mesh_.cellZones()[cellZoneID_];
        const scalarField& V = mesh_.V();
        vectorField& Usource = UEqn.source();
        const vectorField& U = UEqn.psi();
    
        const vector Omega = this->Omega();
    
    andy's avatar
    andy committed
        if (rhs)
    
    andy's avatar
    andy committed
            forAll(cells, i)
            {
                label celli = cells[i];
                Usource[celli] += V[celli]*(Omega ^ U[celli]);
            }
        }
        else
        {
            forAll(cells, i)
            {
                label celli = cells[i];
                Usource[celli] -= V[celli]*(Omega ^ U[celli]);
            }
    
    void Foam::MRFZone::addCoriolis
    (
        const volScalarField& rho,
    
    andy's avatar
    andy committed
        fvVectorMatrix& UEqn,
        const bool rhs
    
    ) const
    {
        if (cellZoneID_ == -1)
        {
            return;
        }
    
        const labelList& cells = mesh_.cellZones()[cellZoneID_];
        const scalarField& V = mesh_.V();
        vectorField& Usource = UEqn.source();
        const vectorField& U = UEqn.psi();
    
        const vector Omega = this->Omega();
    
    andy's avatar
    andy committed
        if (rhs)
    
    andy's avatar
    andy committed
            forAll(cells, i)
            {
                label celli = cells[i];
                Usource[celli] += V[celli]*rho[celli]*(Omega ^ U[celli]);
            }
        }
        else
        {
            forAll(cells, i)
            {
                label celli = cells[i];
                Usource[celli] -= V[celli]*rho[celli]*(Omega ^ U[celli]);
            }
    
    void Foam::MRFZone::makeRelative(volVectorField& U) const
    
        const volVectorField& C = mesh_.C();
    
    
        const vector Omega = this->Omega();
    
    
        const labelList& cells = mesh_.cellZones()[cellZoneID_];
    
        forAll(cells, i)
        {
            label celli = cells[i];
    
    mattijs's avatar
    mattijs committed
            U[celli] -= (Omega ^ (C[celli] - origin_));
    
        volVectorField::Boundary& Ubf = U.boundaryFieldRef();
    
        forAll(includedFaces_, patchi)
        {
            forAll(includedFaces_[patchi], i)
            {
                label patchFacei = includedFaces_[patchi][i];
    
                Ubf[patchi][patchFacei] = Zero;
    
            }
        }
    
        // Excluded patches
        forAll(excludedFaces_, patchi)
        {
            forAll(excludedFaces_[patchi], i)
            {
                label patchFacei = excludedFaces_[patchi][i];
    
    mattijs's avatar
    mattijs committed
                  ^ (C.boundaryField()[patchi][patchFacei] - origin_));
    
    void Foam::MRFZone::makeRelative(surfaceScalarField& phi) const
    {
        makeRelativeRhoFlux(geometricOneField(), phi);
    }
    
    
    void Foam::MRFZone::makeRelative(FieldField<fvsPatchField, scalar>& phi) const
    {
        makeRelativeRhoFlux(oneFieldField(), phi);
    }
    
    
    
    void Foam::MRFZone::makeRelative(Field<scalar>& phi, const label patchi) const
    {
        makeRelativeRhoFlux(oneField(), phi, patchi);
    }
    
    
    
    void Foam::MRFZone::makeRelative
    (
        const surfaceScalarField& rho,
        surfaceScalarField& phi
    ) const
    {
        makeRelativeRhoFlux(rho, phi);
    }
    
    
    
    void Foam::MRFZone::makeAbsolute(volVectorField& U) const
    
        const volVectorField& C = mesh_.C();
    
    
        const vector Omega = this->Omega();
    
    
        const labelList& cells = mesh_.cellZones()[cellZoneID_];
    
        forAll(cells, i)
        {
            label celli = cells[i];
    
    mattijs's avatar
    mattijs committed
            U[celli] += (Omega ^ (C[celli] - origin_));
    
        volVectorField::Boundary& Ubf = U.boundaryFieldRef();
    
        forAll(includedFaces_, patchi)
        {
            forAll(includedFaces_[patchi], i)
            {
                label patchFacei = includedFaces_[patchi][i];
    
    mattijs's avatar
    mattijs committed
                    (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
    
            }
        }
    
        // Excluded patches
        forAll(excludedFaces_, patchi)
        {
            forAll(excludedFaces_[patchi], i)
            {
                label patchFacei = excludedFaces_[patchi][i];
    
    mattijs's avatar
    mattijs committed
                    (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
    
    void Foam::MRFZone::makeAbsolute(surfaceScalarField& phi) const
    
        makeAbsoluteRhoFlux(geometricOneField(), phi);
    
    void Foam::MRFZone::makeAbsolute
    
    (
        const surfaceScalarField& rho,
        surfaceScalarField& phi
    ) const
    {
    
    void Foam::MRFZone::correctBoundaryVelocity(volVectorField& U) const
    {
    
        const vector Omega = this->Omega();
    
    mattijs's avatar
    mattijs committed
        // Included patches
    
        volVectorField::Boundary& Ubf = U.boundaryFieldRef();
    
    mattijs's avatar
    mattijs committed
        forAll(includedFaces_, patchi)
    
    mattijs's avatar
    mattijs committed
            const vectorField& patchC = mesh_.Cf().boundaryField()[patchi];
    
    
            vectorField pfld(Ubf[patchi]);
    
    mattijs's avatar
    mattijs committed
    
            forAll(includedFaces_[patchi], i)
            {
                label patchFacei = includedFaces_[patchi][i];
    
    
    mattijs's avatar
    mattijs committed
                pfld[patchFacei] = (Omega ^ (patchC[patchFacei] - origin_));
    
    andy's avatar
    andy committed
    void Foam::MRFZone::writeData(Ostream& os) const
    
    andy's avatar
    andy committed
        os  << nl;
    
        os.beginBlock(name_);
    
        os.writeEntry("active", active_);
        os.writeEntry("cellZone", cellZoneName_);
        os.writeEntry("origin", origin_);
        os.writeEntry("axis", axis_);
    
    andy's avatar
    andy committed
        omega_->writeData(os);
    
    andy's avatar
    andy committed
        if (excludedPatchNames_.size())
    
            os.writeEntry("nonRotatingPatches", excludedPatchNames_);
    
    andy's avatar
    andy committed
    }
    
    andy's avatar
    andy committed
    
    bool Foam::MRFZone::read(const dictionary& dict)
    {
        coeffs_ = dict;
    
    
        coeffs_.readIfPresent("active", active_);
    
        if (!active_)
        {
            cellZoneID_ = -1;
            return true;
        }
    
        coeffs_.readIfPresent("nonRotatingPatches", excludedPatchNames_);
    
        origin_ = coeffs_.get<vector>("origin");
        axis_ = coeffs_.get<vector>("axis").normalise();
        omega_.reset(Function1<scalar>::New("omega", coeffs_, &mesh_));
    
        const word oldCellZoneName = cellZoneName_;
        if (cellZoneName_ == word::null)
        {
            coeffs_.readEntry("cellZone", cellZoneName_);
        }
        else
        {
            coeffs_.readIfPresent("cellZone", cellZoneName_);
        }
    
        if (cellZoneID_ == -1 || oldCellZoneName != cellZoneName_)
        {
            cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
    
            const labelHashSet excludedPatchSet
            (
                mesh_.boundaryMesh().patchSet(excludedPatchNames_)
            );
    
            excludedPatchLabels_.setSize(excludedPatchSet.size());
    
            label i = 0;
            for (const label patchi : excludedPatchSet)
            {
                excludedPatchLabels_[i++] = patchi;
            }
    
            bool cellZoneFound = (cellZoneID_ != -1);
    
            reduce(cellZoneFound, orOp<bool>());
    
            if (!cellZoneFound)
            {
                FatalErrorInFunction
                    << "cannot find MRF cellZone " << cellZoneName_
                    << exit(FatalError);
            }
    
            setMRFFaces();
        }
    
    andy's avatar
    andy committed
    
        return true;
    
    void Foam::MRFZone::update()
    {
        if (mesh_.topoChanging())
        {
            setMRFFaces();
        }
    }
    
    
    
    // ************************************************************************* //