initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / interpolation / surfaceInterpolation / limitedSchemes / blended / blended.H
blob04521df64d9145e7866a9a0f82759e1dd4bbf705
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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::blended
28 Description
29     linear/upwind blended differencing scheme.
31 SourceFiles
32     blended.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef blended_H
37 #define blended_H
39 #include "limitedSurfaceInterpolationScheme.H"
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 namespace Foam
46 /*---------------------------------------------------------------------------*\
47                            Class blended Declaration
48 \*---------------------------------------------------------------------------*/
50 template<class Type>
51 class blended
53     public limitedSurfaceInterpolationScheme<Type>
55     // Private data
57         const scalar blendingFactor_;
60     // Private Member Functions
62         //- Disallow default bitwise copy construct
63         blended(const blended&);
65         //- Disallow default bitwise assignment
66         void operator=(const blended&);
69 public:
71     //- Runtime type information
72     TypeName("blended");
75     // Constructors
77         //- Construct from mesh, faceFlux and blendingFactor
78         blended
79         (
80             const fvMesh& mesh,
81             const surfaceScalarField& faceFlux,
82             const scalar blendingFactor
83         )
84         :
85             limitedSurfaceInterpolationScheme<Type>(mesh, faceFlux),
86             blendingFactor_(blendingFactor)
87         {}
89         //- Construct from mesh and Istream. 
90         //  The name of the flux field is read from the Istream and looked-up
91         //  from the mesh objectRegistry
92         blended
93         (
94             const fvMesh& mesh,
95             Istream& is
96         )
97         :
98             limitedSurfaceInterpolationScheme<Type>(mesh, is),
99             blendingFactor_(readScalar(is))
100         {}
102         //- Construct from mesh, faceFlux and Istream
103         blended
104         (
105             const fvMesh& mesh,
106             const surfaceScalarField& faceFlux,
107             Istream& is
108         )
109         :
110             limitedSurfaceInterpolationScheme<Type>(mesh, faceFlux),
111             blendingFactor_(readScalar(is))
112         {}
115     // Member Functions
117         //- Return the interpolation limiter
118         virtual tmp<surfaceScalarField> limiter
119         (
120             const GeometricField<Type, fvPatchField, volMesh>&
121         ) const
122         {
123             return tmp<surfaceScalarField>
124             (
125                 new surfaceScalarField
126                 (
127                     IOobject
128                     (
129                         "blendedLimiter",
130                         this->mesh().time().timeName(),
131                         this->mesh()
132                     ),
133                     this->mesh(),
134                     dimensionedScalar
135                     (
136                         "blendedLimiter",
137                         dimless,
138                         1 - blendingFactor_
139                     )
140                 )
141             );
142         }
144         //- Return the interpolation weighting factors
145         virtual tmp<surfaceScalarField> weights
146         (
147             const GeometricField<Type, fvPatchField, volMesh>& vf
148         ) const
149         {
150             return
151                 blendingFactor_*this->mesh().surfaceInterpolation::weights()
152               + (1 - blendingFactor_)*pos(this->faceFlux_);
153         }
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 } // End namespace Foam
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 #endif
165 // ************************************************************************* //