initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / interpolation / surfaceInterpolation / schemes / localBlended / localBlended.H
blobd4dcfd4e18021c3a49004759c83273f848ff1546
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::localBlended
28 Description
29     Two-scheme localBlended differencing scheme.
31 SourceFiles
32     localBlended.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef localBlended_H
37 #define localBlended_H
39 #include "surfaceInterpolationScheme.H"
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 namespace Foam
46 /*---------------------------------------------------------------------------*\
47                            Class localBlended Declaration
48 \*---------------------------------------------------------------------------*/
50 template<class Type>
51 class localBlended
53     public surfaceInterpolationScheme<Type>
55     // Private Member Functions
57         //- Scheme 1
58         tmp<surfaceInterpolationScheme<Type> > tScheme1_;
60         //- Scheme 2
61         tmp<surfaceInterpolationScheme<Type> > tScheme2_;
64         //- Disallow default bitwise copy construct
65         localBlended(const localBlended&);
67         //- Disallow default bitwise assignment
68         void operator=(const localBlended&);
71 public:
73     //- Runtime type information
74     TypeName("localBlended");
77     // Constructors
79         //- Construct from mesh and Istream. 
80         //  The name of the flux field is read from the Istream and looked-up
81         //  from the mesh objectRegistry
82         localBlended
83         (
84             const fvMesh& mesh,
85             Istream& is
86         )
87         :
88             surfaceInterpolationScheme<Type>(mesh),
89             tScheme1_
90             (
91                 surfaceInterpolationScheme<Type>::New(mesh, is)
92             ),
93             tScheme2_
94             (
95                 surfaceInterpolationScheme<Type>::New(mesh, is)
96             )
97         {}
99         //- Construct from mesh, faceFlux and Istream
100         localBlended
101         (
102             const fvMesh& mesh,
103             const surfaceScalarField& faceFlux,
104             Istream& is
105         )
106         :
107             surfaceInterpolationScheme<Type>(mesh),
108             tScheme1_
109             (
110                 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
111             ),
112             tScheme2_
113             (
114                 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
115             )
116         {}
119     // Member Functions
121         //- Return the interpolation weighting factors
122         tmp<surfaceScalarField> weights
123         (
124             const GeometricField<Type, fvPatchField, volMesh>& vf
125         ) const
126         {
127             const surfaceScalarField& blendingFactor =
128                 this->mesh().objectRegistry::
129                 lookupObject<const surfaceScalarField>
130                 (
131                     word(vf.name() + "BlendingFactor")
132                 );
134             return
135                 blendingFactor*tScheme1_().weights(vf)
136               + (scalar(1) - blendingFactor)*tScheme2_().weights(vf);
137         }
139         //- Return the face-interpolate of the given cell field
140         //  with explicit correction
141         tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
142         interpolate(const GeometricField<Type, fvPatchField, volMesh>& vf) const
143         {
144             const surfaceScalarField& blendingFactor = 
145             (
146                 this->mesh().objectRegistry::
147                 lookupObject<const surfaceScalarField> 
148                 (
149                     word(vf.name() + "BlendingFactor")
150                 )
151             );
153             return
154                 blendingFactor*tScheme1_().interpolate(vf)
155               + (scalar(1) - blendingFactor)*tScheme2_().interpolate(vf);
156         }
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 } // End namespace Foam
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 #endif
168 // ************************************************************************* //