FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

outletStabilised.H

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Class
00025     Foam::outletStabilised
00026 
00027 Description
00028     Outlet-stabilised interpolation scheme which applies upwind differencing
00029     to the faces of the cells adjacent to outlets.
00030 
00031     This is particularly useful to stabilise the velocity at entrainment
00032     boundaries for LES cases using linear or other centred differencing
00033     schemes.
00034 
00035 SourceFiles
00036     outletStabilised.C
00037 
00038 \*---------------------------------------------------------------------------*/
00039 
00040 #ifndef outletStabilised_H
00041 #define outletStabilised_H
00042 
00043 #include <finiteVolume/surfaceInterpolationScheme.H>
00044 #include <finiteVolume/skewCorrectionVectors.H>
00045 #include <finiteVolume/linear.H>
00046 #include <finiteVolume/gaussGrad.H>
00047 #include <finiteVolume/mixedFvPatchField.H>
00048 #include <finiteVolume/directionMixedFvPatchField.H>
00049 
00050 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00051 
00052 namespace Foam
00053 {
00054 
00055 /*---------------------------------------------------------------------------*\
00056                            Class outletStabilised Declaration
00057 \*---------------------------------------------------------------------------*/
00058 
00059 template<class Type>
00060 class outletStabilised
00061 :
00062     public surfaceInterpolationScheme<Type>
00063 {
00064     // Private member data
00065 
00066         const surfaceScalarField& faceFlux_;
00067         tmp<surfaceInterpolationScheme<Type> > tScheme_;
00068 
00069 
00070     // Private Member Functions
00071 
00072         //- Disallow default bitwise copy construct
00073         outletStabilised(const outletStabilised&);
00074 
00075         //- Disallow default bitwise assignment
00076         void operator=(const outletStabilised&);
00077 
00078 
00079 public:
00080 
00081     //- Runtime type information
00082     TypeName("outletStabilised");
00083 
00084 
00085     // Constructors
00086 
00087         //- Construct from mesh and Istream
00088         outletStabilised
00089         (
00090             const fvMesh& mesh,
00091             Istream& is
00092         )
00093         :
00094             surfaceInterpolationScheme<Type>(mesh),
00095             faceFlux_
00096             (
00097                 mesh.lookupObject<surfaceScalarField>
00098                 (
00099                     word(is)
00100                 )
00101             ),
00102             tScheme_
00103             (
00104                 surfaceInterpolationScheme<Type>::New(mesh, faceFlux_, is)
00105             )
00106         {}
00107 
00108 
00109         //- Construct from mesh, faceFlux and Istream
00110         outletStabilised
00111         (
00112             const fvMesh& mesh,
00113             const surfaceScalarField& faceFlux,
00114             Istream& is
00115         )
00116         :
00117             surfaceInterpolationScheme<Type>(mesh),
00118             faceFlux_(faceFlux),
00119             tScheme_
00120             (
00121                 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
00122             )
00123         {}
00124 
00125 
00126     // Member Functions
00127 
00128         //- Return the interpolation weighting factors
00129         tmp<surfaceScalarField> weights
00130         (
00131             const GeometricField<Type, fvPatchField, volMesh>& vf
00132         ) const
00133         {
00134             tmp<surfaceScalarField> tw = tScheme_().weights(vf);
00135             surfaceScalarField& w = tw();
00136 
00137             const fvMesh& mesh_ = this->mesh();
00138             const cellList& cells = mesh_.cells();
00139 
00140             forAll(vf.boundaryField(), patchi)
00141             {
00142                 if
00143                 (
00144                     isA<zeroGradientFvPatchField<Type> >
00145                         (vf.boundaryField()[patchi])
00146                  || isA<mixedFvPatchField<Type> >(vf.boundaryField()[patchi])
00147                  || isA<directionMixedFvPatchField<Type> >
00148                     (vf.boundaryField()[patchi])
00149                 )
00150                 {
00151                     const labelList& pFaceCells =
00152                         mesh_.boundary()[patchi].faceCells();
00153 
00154                     forAll(pFaceCells, pFacei)
00155                     {
00156                         const cell& pFaceCell = cells[pFaceCells[pFacei]];
00157 
00158                         forAll(pFaceCell, fi)
00159                         {
00160                             label facei = pFaceCell[fi];
00161 
00162                             if (mesh_.isInternalFace(facei))
00163                             {
00164                                 // Apply upwind differencing
00165                                 w[facei] = pos(faceFlux_[facei]);
00166                             }
00167                         }
00168                     }
00169                 }
00170             }
00171 
00172             return tw;
00173         }
00174 
00175         //- Return true if this scheme uses an explicit correction
00176         virtual bool corrected() const
00177         {
00178             return tScheme_().corrected();
00179         }
00180 
00181         //- Return the explicit correction to the face-interpolate
00182         //  set to zero on the near-boundary faces where upwinf is applied
00183         virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00184         correction
00185         (
00186             const GeometricField<Type, fvPatchField, volMesh>& vf
00187         ) const
00188         {
00189             if (tScheme_().corrected())
00190             {
00191                 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tcorr =
00192                     tScheme_().correction(vf);
00193 
00194                 GeometricField<Type, fvsPatchField, surfaceMesh>& corr = tcorr();
00195 
00196                 const fvMesh& mesh_ = this->mesh();
00197                 const cellList& cells = mesh_.cells();
00198 
00199                 forAll(vf.boundaryField(), patchi)
00200                 {
00201                     if
00202                     (
00203                         isA<zeroGradientFvPatchField<Type> >
00204                             (vf.boundaryField()[patchi])
00205                      || isA<mixedFvPatchField<Type> >
00206                             (vf.boundaryField()[patchi])
00207                     )
00208                     {
00209                         const labelList& pFaceCells =
00210                             mesh_.boundary()[patchi].faceCells();
00211 
00212                         forAll(pFaceCells, pFacei)
00213                         {
00214                             const cell& pFaceCell = cells[pFaceCells[pFacei]];
00215 
00216                             forAll(pFaceCell, fi)
00217                             {
00218                                 label facei = pFaceCell[fi];
00219 
00220                                 if (mesh_.isInternalFace(facei))
00221                                 {
00222                                     // Remove correction
00223                                     corr[facei] = pTraits<Type>::zero;
00224                                 }
00225                             }
00226                         }
00227                     }
00228                 }
00229 
00230                 return tcorr;
00231             }
00232             else
00233             {
00234                 return
00235                     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >(NULL);
00236             }
00237         }
00238 };
00239 
00240 
00241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00242 
00243 } // End namespace Foam
00244 
00245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00246 
00247 #endif
00248 
00249 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines