initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / interpolation / surfaceInterpolation / schemes / outletStabilised / outletStabilised.H
blobdbe27c59966d11d4b63f59981c0af0d7fb8ba111
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::outletStabilised
28 Description
29     Outlet-stabilised interpolation scheme which applies upwind differencing
30     to the faces of the cells adjacent to outlets.
32     This is particularly useful to stabilise the velocity at entrainment
33     boundaries for LES cases using linear or other centred differencing
34     schemes.
36 SourceFiles
37     outletStabilised.C
39 \*---------------------------------------------------------------------------*/
41 #ifndef outletStabilised_H
42 #define outletStabilised_H
44 #include "surfaceInterpolationScheme.H"
45 #include "skewCorrectionVectors.H"
46 #include "linear.H"
47 #include "gaussGrad.H"
48 #include "mixedFvPatchField.H"
49 #include "directionMixedFvPatchField.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 namespace Foam
56 /*---------------------------------------------------------------------------*\
57                            Class outletStabilised Declaration
58 \*---------------------------------------------------------------------------*/
60 template<class Type>
61 class outletStabilised
63     public surfaceInterpolationScheme<Type>
65     // Private member data
67         const surfaceScalarField& faceFlux_;
68         tmp<surfaceInterpolationScheme<Type> > tScheme_;
71     // Private Member Functions
73         //- Disallow default bitwise copy construct
74         outletStabilised(const outletStabilised&);
76         //- Disallow default bitwise assignment
77         void operator=(const outletStabilised&);
80 public:
82     //- Runtime type information
83     TypeName("outletStabilised");
86     // Constructors
88         //- Construct from mesh and Istream
89         outletStabilised
90         (
91             const fvMesh& mesh,
92             Istream& is
93         )
94         :
95             surfaceInterpolationScheme<Type>(mesh),
96             faceFlux_
97             (
98                 mesh.lookupObject<surfaceScalarField>
99                 (
100                     word(is)
101                 )
102             ),
103             tScheme_
104             (
105                 surfaceInterpolationScheme<Type>::New(mesh, faceFlux_, is)
106             )
107         {}
110         //- Construct from mesh, faceFlux and Istream
111         outletStabilised
112         (
113             const fvMesh& mesh,
114             const surfaceScalarField& faceFlux,
115             Istream& is
116         )
117         :
118             surfaceInterpolationScheme<Type>(mesh),
119             faceFlux_(faceFlux),
120             tScheme_
121             (
122                 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
123             )
124         {}
127     // Member Functions
129         //- Return the interpolation weighting factors
130         tmp<surfaceScalarField> weights
131         (
132             const GeometricField<Type, fvPatchField, volMesh>& vf
133         ) const
134         {
135             tmp<surfaceScalarField> tw = tScheme_().weights(vf);
136             surfaceScalarField& w = tw();
138             const fvMesh& mesh_ = this->mesh();
139             const cellList& cells = mesh_.cells();
141             forAll(vf.boundaryField(), patchi)
142             {
143                 if
144                 (
145                     isA<zeroGradientFvPatchField<Type> >
146                         (vf.boundaryField()[patchi])
147                  || isA<mixedFvPatchField<Type> >(vf.boundaryField()[patchi])
148                  || isA<directionMixedFvPatchField<Type> >
149                     (vf.boundaryField()[patchi])
150                 )
151                 {
152                     const labelList& pFaceCells =
153                         mesh_.boundary()[patchi].faceCells();
155                     forAll(pFaceCells, pFacei)
156                     {
157                         const cell& pFaceCell = cells[pFaceCells[pFacei]];
159                         forAll(pFaceCell, fi)
160                         {
161                             label facei = pFaceCell[fi];
163                             if (mesh_.isInternalFace(facei))
164                             {
165                                 // Apply upwind differencing
166                                 w[facei] = pos(faceFlux_[facei]);
167                             }
168                         }
169                     }
170                 }
171             }
173             return tw;
174         }
176         //- Return true if this scheme uses an explicit correction
177         virtual bool corrected() const
178         {
179             return tScheme_().corrected();
180         }
182         //- Return the explicit correction to the face-interpolate
183         //  set to zero on the near-boundary faces where upwinf is applied
184         virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
185         correction
186         (
187             const GeometricField<Type, fvPatchField, volMesh>& vf
188         ) const
189         {
190             if (tScheme_().corrected())
191             {
192                 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tcorr =
193                     tScheme_().correction(vf);
195                 GeometricField<Type, fvsPatchField, surfaceMesh>& corr = tcorr();
197                 const fvMesh& mesh_ = this->mesh();
198                 const cellList& cells = mesh_.cells();
200                 forAll(vf.boundaryField(), patchi)
201                 {
202                     if
203                     (
204                         isA<zeroGradientFvPatchField<Type> >
205                             (vf.boundaryField()[patchi])
206                      || isA<mixedFvPatchField<Type> >
207                             (vf.boundaryField()[patchi])
208                     )
209                     {
210                         const labelList& pFaceCells =
211                             mesh_.boundary()[patchi].faceCells();
213                         forAll(pFaceCells, pFacei)
214                         {
215                             const cell& pFaceCell = cells[pFaceCells[pFacei]];
217                             forAll(pFaceCell, fi)
218                             {
219                                 label facei = pFaceCell[fi];
221                                 if (mesh_.isInternalFace(facei))
222                                 {
223                                     // Remove correction
224                                     corr[facei] = pTraits<Type>::zero;
225                                 }
226                             }
227                         }
228                     }
229                 }
231                 return tcorr;
232             }
233             else
234             {
235                 return
236                     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >(NULL);
237             }
238         }
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 } // End namespace Foam
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 #endif
250 // ************************************************************************* //