initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / interpolation / surfaceInterpolation / surfaceInterpolation / surfaceInterpolate.C
blob51a1a495146d16d8ce4454f7f27af273e65591c3
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 \*---------------------------------------------------------------------------*/
27 #include "surfaceInterpolate.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 namespace fvc
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 // Return weighting factors for scheme given by name in dictionary
42 template<class Type>
43 tmp<surfaceInterpolationScheme<Type> > scheme
45     const surfaceScalarField& faceFlux,
46     Istream& streamData
49     return surfaceInterpolationScheme<Type>::New
50     (
51         faceFlux.mesh(),
52         faceFlux,
53         streamData
54     );
58 // Return weighting factors for scheme given by name in dictionary
59 template<class Type>
60 tmp<surfaceInterpolationScheme<Type> > scheme
62     const surfaceScalarField& faceFlux,
63     const word& name
66     return surfaceInterpolationScheme<Type>::New
67     (
68         faceFlux.mesh(),
69         faceFlux,
70         faceFlux.mesh().interpolationScheme(name)
71     );
75 // Return weighting factors for scheme given by name in dictionary
76 template<class Type>
77 tmp<surfaceInterpolationScheme<Type> > scheme
79     const fvMesh& mesh,
80     Istream& streamData
83     return surfaceInterpolationScheme<Type>::New
84     (
85         mesh,
86         streamData
87     );
91 // Return weighting factors for scheme given by name in dictionary
92 template<class Type>
93 tmp<surfaceInterpolationScheme<Type> > scheme
95     const fvMesh& mesh,
96     const word& name
99     return surfaceInterpolationScheme<Type>::New
100     (
101         mesh,
102         mesh.interpolationScheme(name)
103     );
107 // Interpolate field onto faces using scheme given by name in dictionary
108 template<class Type>
109 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
110 interpolate
112     const GeometricField<Type, fvPatchField, volMesh>& vf,
113     const surfaceScalarField& faceFlux,
114     Istream& schemeData
117     if (surfaceInterpolation::debug)
118     {
119         Info<< "interpolate"
120             << "(const GeometricField<Type, fvPatchField, volMesh>&, "
121             << "const surfaceScalarField&, Istream&) : "
122             << "interpolating GeometricField<Type, fvPatchField, volMesh> "
123             << endl;
124     }
126     return scheme<Type>(faceFlux, schemeData)().interpolate(vf);
130 // Interpolate field onto faces using scheme given by name in dictionary
131 template<class Type>
132 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
133 interpolate
135     const GeometricField<Type, fvPatchField, volMesh>& vf,
136     const surfaceScalarField& faceFlux,
137     const word& name
140     if (surfaceInterpolation::debug)
141     {
142         Info<< "interpolate"
143             << "(const GeometricField<Type, fvPatchField, volMesh>&, "
144             << "const surfaceScalarField&, const word&) : "
145             << "interpolating GeometricField<Type, fvPatchField, volMesh> "
146             << "using " << name
147             << endl;
148     }
150     return scheme<Type>(faceFlux, name)().interpolate(vf);
153 // Interpolate field onto faces using scheme given by name in dictionary
154 template<class Type>
155 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
156 interpolate
158     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
159     const surfaceScalarField& faceFlux,
160     const word& name
163     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
164         interpolate(tvf(), faceFlux, name);
166     tvf.clear();
168     return tsf;
171 // Interpolate field onto faces using scheme given by name in dictionary
172 template<class Type>
173 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
174 interpolate
176     const GeometricField<Type, fvPatchField, volMesh>& vf,
177     const tmp<surfaceScalarField>& tFaceFlux,
178     const word& name
181     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
182         interpolate(vf, tFaceFlux(), name);
184     tFaceFlux.clear();
186     return tsf;
189 // Interpolate field onto faces using scheme given by name in dictionary
190 template<class Type>
191 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
192 interpolate
194     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
195     const tmp<surfaceScalarField>& tFaceFlux,
196     const word& name
199     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
200         interpolate(tvf(), tFaceFlux(), name);
202     tvf.clear();
203     tFaceFlux.clear();
205     return tsf;
209 // Interpolate field onto faces using scheme given by name in dictionary
210 template<class Type>
211 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
212 interpolate
214     const GeometricField<Type, fvPatchField, volMesh>& vf,
215     Istream& schemeData
218     if (surfaceInterpolation::debug)
219     {
220         Info<< "interpolate"
221             << "(const GeometricField<Type, fvPatchField, volMesh>&, "
222             << "Istream&) : "
223             << "interpolating GeometricField<Type, fvPatchField, volMesh> "
224             << endl;
225     }
227     return scheme<Type>(vf.mesh(), schemeData)().interpolate(vf);
230 // Interpolate field onto faces using scheme given by name in dictionary
231 template<class Type>
232 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
233 interpolate
235     const GeometricField<Type, fvPatchField, volMesh>& vf,
236     const word& name
239     if (surfaceInterpolation::debug)
240     {
241         Info<< "interpolate"
242             << "(const GeometricField<Type, fvPatchField, volMesh>&, "
243             << "const word&) : "
244             << "interpolating GeometricField<Type, fvPatchField, volMesh> "
245             << "using " << name
246             << endl;
247     }
249     return scheme<Type>(vf.mesh(), name)().interpolate(vf);
252 // Interpolate field onto faces using scheme given by name in dictionary
253 template<class Type>
254 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
255 interpolate
257     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
258     const word& name
261     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
262         interpolate(tvf(), name);
264     tvf.clear();
266     return tsf;
270 // Interpolate field onto faces using central differencing
271 template<class Type>
272 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
273 interpolate
275     const GeometricField<Type, fvPatchField, volMesh>& vf
278     if (surfaceInterpolation::debug)
279     {
280         Info<< "interpolate"
281             << "(const GeometricField<Type, fvPatchField, volMesh>&) : "
282             << "interpolating GeometricField<Type, fvPatchField, volMesh> "
283             << "using run-time selected scheme"
284             << endl;
285     }
287     return interpolate(vf, "interpolate(" + vf.name() + ')');
291 // Interpolate field onto faces using central differencing
292 template<class Type>
293 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > 
294 interpolate
296     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
299     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
300         interpolate(tvf());
301     tvf.clear();
302     return tsf;
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 } // End namespace fvc
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312 } // End namespace Foam
314 // ************************************************************************* //