This file is indexed.

/usr/include/freefoam/finiteVolume/outletStabilised.H is in libfreefoam-dev 0.1.0+dfsg-1build1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 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
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 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
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
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
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
     \\/     M anipulation  |
-------------------------------------------------------------------------------
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/>.

Class
    Foam::outletStabilised

Description
    Outlet-stabilised interpolation scheme which applies upwind differencing
    to the faces of the cells adjacent to outlets.

    This is particularly useful to stabilise the velocity at entrainment
    boundaries for LES cases using linear or other centred differencing
    schemes.

SourceFiles
    outletStabilised.C

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

#ifndef outletStabilised_H
#define outletStabilised_H

#include <finiteVolume/surfaceInterpolationScheme.H>
#include <finiteVolume/skewCorrectionVectors.H>
#include <finiteVolume/linear.H>
#include <finiteVolume/gaussGrad.H>
#include <finiteVolume/mixedFvPatchField.H>
#include <finiteVolume/directionMixedFvPatchField.H>

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

namespace Foam
{

/*---------------------------------------------------------------------------*\
                           Class outletStabilised Declaration
\*---------------------------------------------------------------------------*/

template<class Type>
class outletStabilised
:
    public surfaceInterpolationScheme<Type>
{
    // Private member data

        const surfaceScalarField& faceFlux_;
        tmp<surfaceInterpolationScheme<Type> > tScheme_;


    // Private Member Functions

        //- Disallow default bitwise copy construct
        outletStabilised(const outletStabilised&);

        //- Disallow default bitwise assignment
        void operator=(const outletStabilised&);


public:

    //- Runtime type information
    TypeName("outletStabilised");


    // Constructors

        //- Construct from mesh and Istream
        outletStabilised
        (
            const fvMesh& mesh,
            Istream& is
        )
        :
            surfaceInterpolationScheme<Type>(mesh),
            faceFlux_
            (
                mesh.lookupObject<surfaceScalarField>
                (
                    word(is)
                )
            ),
            tScheme_
            (
                surfaceInterpolationScheme<Type>::New(mesh, faceFlux_, is)
            )
        {}


        //- Construct from mesh, faceFlux and Istream
        outletStabilised
        (
            const fvMesh& mesh,
            const surfaceScalarField& faceFlux,
            Istream& is
        )
        :
            surfaceInterpolationScheme<Type>(mesh),
            faceFlux_(faceFlux),
            tScheme_
            (
                surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
            )
        {}


    // Member Functions

        //- Return the interpolation weighting factors
        tmp<surfaceScalarField> weights
        (
            const GeometricField<Type, fvPatchField, volMesh>& vf
        ) const
        {
            tmp<surfaceScalarField> tw = tScheme_().weights(vf);
            surfaceScalarField& w = tw();

            const fvMesh& mesh_ = this->mesh();
            const cellList& cells = mesh_.cells();

            forAll(vf.boundaryField(), patchi)
            {
                if
                (
                    isA<zeroGradientFvPatchField<Type> >
                        (vf.boundaryField()[patchi])
                 || isA<mixedFvPatchField<Type> >(vf.boundaryField()[patchi])
                 || isA<directionMixedFvPatchField<Type> >
                    (vf.boundaryField()[patchi])
                )
                {
                    const labelList& pFaceCells =
                        mesh_.boundary()[patchi].faceCells();

                    forAll(pFaceCells, pFacei)
                    {
                        const cell& pFaceCell = cells[pFaceCells[pFacei]];

                        forAll(pFaceCell, fi)
                        {
                            label facei = pFaceCell[fi];

                            if (mesh_.isInternalFace(facei))
                            {
                                // Apply upwind differencing
                                w[facei] = pos(faceFlux_[facei]);
                            }
                        }
                    }
                }
            }

            return tw;
        }

        //- Return true if this scheme uses an explicit correction
        virtual bool corrected() const
        {
            return tScheme_().corrected();
        }

        //- Return the explicit correction to the face-interpolate
        //  set to zero on the near-boundary faces where upwinf is applied
        virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
        correction
        (
            const GeometricField<Type, fvPatchField, volMesh>& vf
        ) const
        {
            if (tScheme_().corrected())
            {
                tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tcorr =
                    tScheme_().correction(vf);

                GeometricField<Type, fvsPatchField, surfaceMesh>& corr = tcorr();

                const fvMesh& mesh_ = this->mesh();
                const cellList& cells = mesh_.cells();

                forAll(vf.boundaryField(), patchi)
                {
                    if
                    (
                        isA<zeroGradientFvPatchField<Type> >
                            (vf.boundaryField()[patchi])
                     || isA<mixedFvPatchField<Type> >
                            (vf.boundaryField()[patchi])
                    )
                    {
                        const labelList& pFaceCells =
                            mesh_.boundary()[patchi].faceCells();

                        forAll(pFaceCells, pFacei)
                        {
                            const cell& pFaceCell = cells[pFaceCells[pFacei]];

                            forAll(pFaceCell, fi)
                            {
                                label facei = pFaceCell[fi];

                                if (mesh_.isInternalFace(facei))
                                {
                                    // Remove correction
                                    corr[facei] = pTraits<Type>::zero;
                                }
                            }
                        }
                    }
                }

                return tcorr;
            }
            else
            {
                return
                    tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >(NULL);
            }
        }
};


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

} // End namespace Foam

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

#endif

// ************************ vim: set sw=4 sts=4 et: ************************ //