initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / interpolation / surfaceInterpolation / surfaceInterpolationScheme / surfaceInterpolationScheme.C
blob4702ccb5aa7ed7f23c9f0366b658b8269bd6a066
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 Description
26     Abstract base class for surface interpolation schemes.
28 \*---------------------------------------------------------------------------*/
30 #include "surfaceInterpolationScheme.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "coupledFvPatchField.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace Foam
40 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
42 // Return weighting factors for scheme given by name in dictionary
43 template<class Type>
44 tmp<surfaceInterpolationScheme<Type> > surfaceInterpolationScheme<Type>::New
46     const fvMesh& mesh,
47     Istream& schemeData
50     if (schemeData.eof())
51     {
52         FatalIOErrorIn
53         (
54             "surfaceInterpolationScheme<Type>::New(const fvMesh&, Istream&)",
55             schemeData
56         )   << "Discretisation scheme not specified"
57             << endl << endl
58             << "Valid schemes are :" << endl
59             << MeshConstructorTablePtr_->toc()
60             << exit(FatalIOError);
61     }
63     word schemeName(schemeData);
65     if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
66     {
67         Info<< "surfaceInterpolationScheme<Type>::New"
68                "(const fvMesh&, Istream&)"
69                " : discretisation scheme = "
70             << schemeName
71             << endl;
72     }
74     typename MeshConstructorTable::iterator constructorIter =
75         MeshConstructorTablePtr_->find(schemeName);
77     if (constructorIter == MeshConstructorTablePtr_->end())
78     {
79         FatalIOErrorIn
80         (
81             "surfaceInterpolationScheme<Type>::New(const fvMesh&, Istream&)",
82             schemeData
83         )   << "Unknown discretisation scheme " << schemeName
84             << endl << endl
85             << "Valid schemes are :" << endl
86             << MeshConstructorTablePtr_->toc()
87             << exit(FatalIOError);
88     }
90     return constructorIter()(mesh, schemeData);
94 // Return weighting factors for scheme given by name in dictionary
95 template<class Type>
96 tmp<surfaceInterpolationScheme<Type> > surfaceInterpolationScheme<Type>::New
98     const fvMesh& mesh,
99     const surfaceScalarField& faceFlux,
100     Istream& schemeData
103     if (schemeData.eof())
104     {
105         FatalIOErrorIn
106         (
107             "surfaceInterpolationScheme<Type>::New"
108             "(const fvMesh&, const surfaceScalarField&, Istream&)",
109             schemeData
110         )   << "Discretisation scheme not specified"
111             << endl << endl
112             << "Valid schemes are :" << endl
113             << MeshConstructorTablePtr_->toc()
114             << exit(FatalIOError);
115     }
117     word schemeName(schemeData);
119     if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
120     {
121         Info<< "surfaceInterpolationScheme<Type>::New"
122                "(const fvMesh&, const surfaceScalarField&, Istream&)"
123                " : discretisation scheme = "
124             << schemeName
125             << endl;
126     }
128     typename MeshFluxConstructorTable::iterator constructorIter =
129         MeshFluxConstructorTablePtr_->find(schemeName);
131     if (constructorIter == MeshFluxConstructorTablePtr_->end())
132     {
133         FatalIOErrorIn
134         (
135             "surfaceInterpolationScheme<Type>::New"
136             "(const fvMesh&, const surfaceScalarField&, Istream&)",
137             schemeData
138         )   << "Unknown discretisation scheme " << schemeName
139             << endl << endl
140             << "Valid schemes are :" << endl
141             << MeshFluxConstructorTablePtr_->toc()
142             << exit(FatalIOError);
143     }
145     return constructorIter()(mesh, faceFlux, schemeData);
149 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
151 template<class Type>
152 surfaceInterpolationScheme<Type>::~surfaceInterpolationScheme()
156 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
158 //- Return the face-interpolate of the given cell field
159 //  with the given owner and neighbour weighting factors
160 template<class Type>
161 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
162 surfaceInterpolationScheme<Type>::interpolate
164     const GeometricField<Type, fvPatchField, volMesh>& vf,
165     const tmp<surfaceScalarField>& tlambdas,
166     const tmp<surfaceScalarField>& tys
169     if (surfaceInterpolation::debug)
170     {
171         Info<< "surfaceInterpolationScheme<Type>::uncorrectedInterpolate"
172                "(const GeometricField<Type, fvPatchField, volMesh>&, "
173                "const tmp<surfaceScalarField>&, "
174                "const tmp<surfaceScalarField>&) : "
175                "interpolating volTypeField from cells to faces "
176                "without explicit correction"
177             << endl;
178     }
180     const surfaceScalarField& lambdas = tlambdas();
181     const surfaceScalarField& ys = tys();
183     const Field<Type>& vfi = vf.internalField();
184     const scalarField& lambda = lambdas.internalField();
185     const scalarField& y = ys.internalField();
187     const fvMesh& mesh = vf.mesh();
188     const unallocLabelList& P = mesh.owner();
189     const unallocLabelList& N = mesh.neighbour();
191     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf
192     (
193         new GeometricField<Type, fvsPatchField, surfaceMesh>
194         (
195             IOobject
196             (
197                 "interpolate("+vf.name()+')',
198                 vf.instance(),
199                 vf.db()
200             ),
201             mesh,
202             vf.dimensions()
203         )
204     );
205     GeometricField<Type, fvsPatchField, surfaceMesh>& sf = tsf();
207     Field<Type>& sfi = sf.internalField();
209     for (register label fi=0; fi<P.size(); fi++)
210     {
211         sfi[fi] = lambda[fi]*vfi[P[fi]] + y[fi]*vfi[N[fi]];
212     }
215     // Interpolate across coupled patches using given lambdas and ys
217     forAll (lambdas.boundaryField(), pi)
218     {
219         const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
220         const fvsPatchScalarField& pY = ys.boundaryField()[pi];
222         if (vf.boundaryField()[pi].coupled())
223         {
224             sf.boundaryField()[pi] =
225                 pLambda*vf.boundaryField()[pi].patchInternalField()
226               + pY*vf.boundaryField()[pi].patchNeighbourField();
227         }
228         else
229         {
230             sf.boundaryField()[pi] = vf.boundaryField()[pi];
231         }
232     }
234     tlambdas.clear();
235     tys.clear();
237     return tsf;
241 //- Return the face-interpolate of the given cell field
242 //  with the given weigting factors
243 template<class Type>
244 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
245 surfaceInterpolationScheme<Type>::interpolate
247     const GeometricField<Type, fvPatchField, volMesh>& vf,
248     const tmp<surfaceScalarField>& tlambdas
251     if (surfaceInterpolation::debug)
252     {
253         Info<< "surfaceInterpolationScheme<Type>::interpolate"
254                "(const GeometricField<Type, fvPatchField, volMesh>&, "
255                "const tmp<surfaceScalarField>&) : "
256                "interpolating volTypeField from cells to faces "
257                "without explicit correction"
258             << endl;
259     }
261     const surfaceScalarField& lambdas = tlambdas();
263     const Field<Type>& vfi = vf.internalField();
264     const scalarField& lambda = lambdas.internalField();
266     const fvMesh& mesh = vf.mesh();
267     const unallocLabelList& P = mesh.owner();
268     const unallocLabelList& N = mesh.neighbour();
270     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf
271     (
272         new GeometricField<Type, fvsPatchField, surfaceMesh>
273         (
274             IOobject
275             (
276                 "interpolate("+vf.name()+')',
277                 vf.instance(),
278                 vf.db()
279             ),
280             mesh,
281             vf.dimensions()
282         )
283     );
284     GeometricField<Type, fvsPatchField, surfaceMesh>& sf = tsf();
286     Field<Type>& sfi = sf.internalField();
288     for (register label fi=0; fi<P.size(); fi++)
289     {
290         sfi[fi] = lambda[fi]*(vfi[P[fi]] - vfi[N[fi]]) + vfi[N[fi]];
291     }
293     // Interpolate across coupled patches using given lambdas
295     forAll (lambdas.boundaryField(), pi)
296     {
297         const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
299         if (vf.boundaryField()[pi].coupled())
300         {
301             tsf().boundaryField()[pi] =
302                 pLambda*vf.boundaryField()[pi].patchInternalField()
303              + (1.0 - pLambda)*vf.boundaryField()[pi].patchNeighbourField();
304         }
305         else
306         {
307             sf.boundaryField()[pi] = vf.boundaryField()[pi];
308         }
309     }
311     tlambdas.clear();
313     return tsf;
317 //- Return the face-interpolate of the given cell field
318 //  with explicit correction
319 template<class Type>
320 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
321 surfaceInterpolationScheme<Type>::interpolate
323     const GeometricField<Type, fvPatchField, volMesh>& vf
324 ) const
326     if (surfaceInterpolation::debug)
327     {
328         Info<< "surfaceInterpolationScheme<Type>::interpolate"
329                "(const GeometricField<Type, fvPatchField, volMesh>&) : "
330             << "interpolating volTypeField from cells to faces"
331             << endl;
332     }
334     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf
335         = interpolate(vf, weights(vf));
337     if (corrected())
338     {
339         tsf() += correction(vf);
340     }
342     return tsf;
346 //- Return the face-interpolate of the given cell field
347 //  with explicit correction
348 template<class Type>
349 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
350 surfaceInterpolationScheme<Type>::interpolate
352     const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
353 ) const
355     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tinterpVf
356         = interpolate(tvf());
357     tvf.clear();
358     return tinterpVf;
362 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 } // End namespace Foam
366 // ************************************************************************* //