BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / finiteVolume / interpolation / surfaceInterpolation / schemes / localBlended / localBlended.H
blob5f49be98fd2348f75b2586ccf225c602fdaa13b3
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::localBlended
27 Description
28     Two-scheme localBlended differencing scheme.
30 SourceFiles
31     localBlended.C
33 \*---------------------------------------------------------------------------*/
35 #ifndef localBlended_H
36 #define localBlended_H
38 #include "surfaceInterpolationScheme.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 namespace Foam
45 /*---------------------------------------------------------------------------*\
46                            Class localBlended Declaration
47 \*---------------------------------------------------------------------------*/
49 template<class Type>
50 class localBlended
52     public surfaceInterpolationScheme<Type>
54     // Private Member Functions
56         //- Scheme 1
57         tmp<surfaceInterpolationScheme<Type> > tScheme1_;
59         //- Scheme 2
60         tmp<surfaceInterpolationScheme<Type> > tScheme2_;
63         //- Disallow default bitwise copy construct
64         localBlended(const localBlended&);
66         //- Disallow default bitwise assignment
67         void operator=(const localBlended&);
70 public:
72     //- Runtime type information
73     TypeName("localBlended");
76     // Constructors
78         //- Construct from mesh and Istream.
79         //  The name of the flux field is read from the Istream and looked-up
80         //  from the mesh objectRegistry
81         localBlended
82         (
83             const fvMesh& mesh,
84             Istream& is
85         )
86         :
87             surfaceInterpolationScheme<Type>(mesh),
88             tScheme1_
89             (
90                 surfaceInterpolationScheme<Type>::New(mesh, is)
91             ),
92             tScheme2_
93             (
94                 surfaceInterpolationScheme<Type>::New(mesh, is)
95             )
96         {}
98         //- Construct from mesh, faceFlux and Istream
99         localBlended
100         (
101             const fvMesh& mesh,
102             const surfaceScalarField& faceFlux,
103             Istream& is
104         )
105         :
106             surfaceInterpolationScheme<Type>(mesh),
107             tScheme1_
108             (
109                 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
110             ),
111             tScheme2_
112             (
113                 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
114             )
115         {}
118     // Member Functions
120         //- Return the interpolation weighting factors
121         tmp<surfaceScalarField> weights
122         (
123             const GeometricField<Type, fvPatchField, volMesh>& vf
124         ) const
125         {
126             const surfaceScalarField& blendingFactor =
127                 this->mesh().objectRegistry::template
128                 lookupObject<const surfaceScalarField>
129                 (
130                     word(vf.name() + "BlendingFactor")
131                 );
133             return
134                 blendingFactor*tScheme1_().weights(vf)
135               + (scalar(1) - blendingFactor)*tScheme2_().weights(vf);
136         }
138         //- Return the face-interpolate of the given cell field
139         //  with explicit correction
140         tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
141         interpolate(const GeometricField<Type, fvPatchField, volMesh>& vf) const
142         {
143             const surfaceScalarField& blendingFactor =
144             (
145                 this->mesh().objectRegistry::template
146                 lookupObject<const surfaceScalarField>
147                 (
148                     word(vf.name() + "BlendingFactor")
149                 )
150             );
152             return
153                 blendingFactor*tScheme1_().interpolate(vf)
154               + (scalar(1) - blendingFactor)*tScheme2_().interpolate(vf);
155         }
158         //- Return true if this scheme uses an explicit correction
159         virtual bool corrected() const
160         {
161             return tScheme1_().corrected() || tScheme2_().corrected();
162         }
165         //- Return the explicit correction to the face-interpolate
166         //  for the given field
167         virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
168         correction
169         (
170             const GeometricField<Type, fvPatchField, volMesh>& vf
171         ) const
172         {
173             const surfaceScalarField& blendingFactor =
174                 this->mesh().objectRegistry::template
175                 lookupObject<const surfaceScalarField>
176                 (
177                     word(vf.name() + "BlendingFactor")
178                 );
180             if (tScheme1_().corrected())
181             {
182                 if (tScheme2_().corrected())
183                 {
184                     return
185                     (
186                         blendingFactor
187                       * tScheme1_().correction(vf)
188                       + (scalar(1.0) - blendingFactor)
189                       * tScheme2_().correction(vf)
190                     );
191                 }
192                 else
193                 {
194                     return
195                     (
196                         blendingFactor
197                       * tScheme1_().correction(vf)
198                     );
199                 }
200             }
201             else if (tScheme2_().corrected())
202             {
203                 return
204                 (
205                     (scalar(1.0) - blendingFactor)
206                   * tScheme2_().correction(vf)
207                 );
208             }
209             else
210             {
211                 return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
212                 (
213                     NULL
214                 );
215             }
216         }
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 } // End namespace Foam
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 #endif
228 // ************************************************************************* //