1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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
26 Foam::outletStabilised
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
39 \*---------------------------------------------------------------------------*/
41 #ifndef outletStabilised_H
42 #define outletStabilised_H
44 #include "surfaceInterpolationScheme.H"
45 #include "skewCorrectionVectors.H"
47 #include "gaussGrad.H"
48 #include "mixedFvPatchField.H"
49 #include "directionMixedFvPatchField.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 /*---------------------------------------------------------------------------*\
57 Class outletStabilised Declaration
58 \*---------------------------------------------------------------------------*/
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&);
82 //- Runtime type information
83 TypeName("outletStabilised");
88 //- Construct from mesh and Istream
95 surfaceInterpolationScheme<Type>(mesh),
98 mesh.lookupObject<surfaceScalarField>
105 surfaceInterpolationScheme<Type>::New(mesh, faceFlux_, is)
110 //- Construct from mesh, faceFlux and Istream
114 const surfaceScalarField& faceFlux,
118 surfaceInterpolationScheme<Type>(mesh),
122 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
129 //- Return the interpolation weighting factors
130 tmp<surfaceScalarField> weights
132 const GeometricField<Type, fvPatchField, volMesh>& vf
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)
145 isA<zeroGradientFvPatchField<Type> >
146 (vf.boundaryField()[patchi])
147 || isA<mixedFvPatchField<Type> >(vf.boundaryField()[patchi])
148 || isA<directionMixedFvPatchField<Type> >
149 (vf.boundaryField()[patchi])
152 const labelList& pFaceCells =
153 mesh_.boundary()[patchi].faceCells();
155 forAll(pFaceCells, pFacei)
157 const cell& pFaceCell = cells[pFaceCells[pFacei]];
159 forAll(pFaceCell, fi)
161 label facei = pFaceCell[fi];
163 if (mesh_.isInternalFace(facei))
165 // Apply upwind differencing
166 w[facei] = pos(faceFlux_[facei]);
176 //- Return true if this scheme uses an explicit correction
177 virtual bool corrected() const
179 return tScheme_().corrected();
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> >
187 const GeometricField<Type, fvPatchField, volMesh>& vf
190 if (tScheme_().corrected())
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)
204 isA<zeroGradientFvPatchField<Type> >
205 (vf.boundaryField()[patchi])
206 || isA<mixedFvPatchField<Type> >
207 (vf.boundaryField()[patchi])
210 const labelList& pFaceCells =
211 mesh_.boundary()[patchi].faceCells();
213 forAll(pFaceCells, pFacei)
215 const cell& pFaceCell = cells[pFaceCells[pFacei]];
217 forAll(pFaceCell, fi)
219 label facei = pFaceCell[fi];
221 if (mesh_.isInternalFace(facei))
224 corr[facei] = pTraits<Type>::zero;
236 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >(NULL);
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 } // End namespace Foam
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 // ************************************************************************* //