LimitedScheme.C 6.07 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2016 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 "volFields.H"
#include "surfaceFields.H"
#include "fvcGrad.H"
#include "coupledFvPatchFields.H"

31
// * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //
32 33

template<class Type, class Limiter, template<class> class LimitFunc>
34
void Foam::LimitedScheme<Type, Limiter, LimitFunc>::calcLimiter
35
(
36 37
    const GeometricField<Type, fvPatchField, volMesh>& phi,
    surfaceScalarField& limiterField
38 39
) const
{
40 41 42 43 44
    typedef GeometricField<typename Limiter::phiType, fvPatchField, volMesh>
        VolFieldType;

    typedef GeometricField<typename Limiter::gradPhiType, fvPatchField, volMesh>
        GradVolFieldType;
45

46
    const fvMesh& mesh = this->mesh();
47

48 49
    tmp<VolFieldType> tlPhi = LimitFunc<Type>()(phi);
    const VolFieldType& lPhi = tlPhi();
50

51 52
    tmp<GradVolFieldType> tgradc(fvc::grad(lPhi));
    const GradVolFieldType& gradc = tgradc();
53 54 55

    const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights();

56 57
    const labelUList& owner = mesh.owner();
    const labelUList& neighbour = mesh.neighbour();
58 59 60

    const vectorField& C = mesh.C();

61
    scalarField& pLim = limiterField.primitiveFieldRef();
62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79

    forAll(pLim, face)
    {
        label own = owner[face];
        label nei = neighbour[face];

        pLim[face] = Limiter::limiter
        (
            CDweights[face],
            this->faceFlux_[face],
            lPhi[own],
            lPhi[nei],
            gradc[own],
            gradc[nei],
            C[nei] - C[own]
        );
    }

80
    surfaceScalarField::Boundary& bLim = limiterField.boundaryFieldRef();
81 82 83 84 85 86 87 88 89 90

    forAll(bLim, patchi)
    {
        scalarField& pLim = bLim[patchi];

        if (bLim[patchi].coupled())
        {
            const scalarField& pCDweights = CDweights.boundaryField()[patchi];
            const scalarField& pFaceFlux =
                this->faceFlux_.boundaryField()[patchi];
91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107

            const Field<typename Limiter::phiType> plPhiP
            (
                lPhi.boundaryField()[patchi].patchInternalField()
            );
            const Field<typename Limiter::phiType> plPhiN
            (
                lPhi.boundaryField()[patchi].patchNeighbourField()
            );
            const Field<typename Limiter::gradPhiType> pGradcP
            (
                gradc.boundaryField()[patchi].patchInternalField()
            );
            const Field<typename Limiter::gradPhiType> pGradcN
            (
                gradc.boundaryField()[patchi].patchNeighbourField()
            );
108 109

            // Build the d-vectors
Henry's avatar
Henry committed
110
            vectorField pd(CDweights.boundaryField()[patchi].patch().delta());
111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130

            forAll(pLim, face)
            {
                pLim[face] = Limiter::limiter
                (
                    pCDweights[face],
                    pFaceFlux[face],
                    plPhiP[face],
                    plPhiN[face],
                    pGradcP[face],
                    pGradcN[face],
                    pd[face]
                );
            }
        }
        else
        {
            pLim = 1.0;
        }
    }
131 132

    limiterField.setOriented();
133 134 135 136
}


// * * * * * * * * * * * * Public Member Functions  * * * * * * * * * * * * //
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 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202
template<class Type, class Limiter, template<class> class LimitFunc>
Foam::tmp<Foam::surfaceScalarField>
Foam::LimitedScheme<Type, Limiter, LimitFunc>::limiter
(
    const GeometricField<Type, fvPatchField, volMesh>& phi
) const
{
    const fvMesh& mesh = this->mesh();

    const word limiterFieldName(type() + "Limiter(" + phi.name() + ')');

    if (this->mesh().cache("limiter"))
    {
        if (!mesh.foundObject<surfaceScalarField>(limiterFieldName))
        {
            surfaceScalarField* limiterField
            (
                new surfaceScalarField
                (
                    IOobject
                    (
                        limiterFieldName,
                        mesh.time().timeName(),
                        mesh,
                        IOobject::NO_READ,
                        IOobject::NO_WRITE
                    ),
                    mesh,
                    dimless
                )
            );

            mesh.objectRegistry::store(limiterField);
        }

        surfaceScalarField& limiterField =
            const_cast<surfaceScalarField&>
            (
                mesh.lookupObject<surfaceScalarField>
                (
                    limiterFieldName
                )
            );

        calcLimiter(phi, limiterField);

        return limiterField;
    }
    else
    {
        tmp<surfaceScalarField> tlimiterField
        (
            new surfaceScalarField
            (
                IOobject
                (
                    limiterFieldName,
                    mesh.time().timeName(),
                    mesh
                ),
                mesh,
                dimless
            )
        );

203
        calcLimiter(phi, tlimiterField.ref());
204 205 206

        return tlimiterField;
    }
207 208 209 210
}


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