1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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/>.
25 Abstract base class for surface interpolation schemes.
27 \*---------------------------------------------------------------------------*/
29 #include "surfaceInterpolationScheme.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "coupledFvPatchField.H"
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
41 // Return weighting factors for scheme given by name in dictionary
43 tmp<surfaceInterpolationScheme<Type> > surfaceInterpolationScheme<Type>::New
53 "surfaceInterpolationScheme<Type>::New(const fvMesh&, Istream&)",
55 ) << "Discretisation scheme not specified"
57 << "Valid schemes are :" << endl
58 << MeshConstructorTablePtr_->sortedToc()
59 << exit(FatalIOError);
62 const word schemeName(schemeData);
64 if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
66 Info<< "surfaceInterpolationScheme<Type>::New"
67 "(const fvMesh&, Istream&)"
68 " : discretisation scheme = "
73 typename MeshConstructorTable::iterator constructorIter =
74 MeshConstructorTablePtr_->find(schemeName);
76 if (constructorIter == MeshConstructorTablePtr_->end())
80 "surfaceInterpolationScheme<Type>::New(const fvMesh&, Istream&)",
82 ) << "Unknown discretisation scheme "
83 << schemeName << nl << nl
84 << "Valid schemes are :" << endl
85 << MeshConstructorTablePtr_->sortedToc()
86 << exit(FatalIOError);
89 return constructorIter()(mesh, schemeData);
93 // Return weighting factors for scheme given by name in dictionary
95 tmp<surfaceInterpolationScheme<Type> > surfaceInterpolationScheme<Type>::New
98 const surfaceScalarField& faceFlux,
102 if (schemeData.eof())
106 "surfaceInterpolationScheme<Type>::New"
107 "(const fvMesh&, const surfaceScalarField&, Istream&)",
109 ) << "Discretisation scheme not specified"
111 << "Valid schemes are :" << endl
112 << MeshConstructorTablePtr_->sortedToc()
113 << exit(FatalIOError);
116 const word schemeName(schemeData);
118 if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
120 Info<< "surfaceInterpolationScheme<Type>::New"
121 "(const fvMesh&, const surfaceScalarField&, Istream&)"
122 " : discretisation scheme = "
127 typename MeshFluxConstructorTable::iterator constructorIter =
128 MeshFluxConstructorTablePtr_->find(schemeName);
130 if (constructorIter == MeshFluxConstructorTablePtr_->end())
134 "surfaceInterpolationScheme<Type>::New"
135 "(const fvMesh&, const surfaceScalarField&, Istream&)",
137 ) << "Unknown discretisation scheme "
138 << schemeName << nl << nl
139 << "Valid schemes are :" << endl
140 << MeshFluxConstructorTablePtr_->sortedToc()
141 << exit(FatalIOError);
144 return constructorIter()(mesh, faceFlux, schemeData);
148 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
151 surfaceInterpolationScheme<Type>::~surfaceInterpolationScheme()
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157 //- Return the face-interpolate of the given cell field
158 // with the given owner and neighbour weighting factors
160 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
161 surfaceInterpolationScheme<Type>::interpolate
163 const GeometricField<Type, fvPatchField, volMesh>& vf,
164 const tmp<surfaceScalarField>& tlambdas,
165 const tmp<surfaceScalarField>& tys
168 if (surfaceInterpolation::debug)
170 Info<< "surfaceInterpolationScheme<Type>::uncorrectedInterpolate"
171 "(const GeometricField<Type, fvPatchField, volMesh>&, "
172 "const tmp<surfaceScalarField>&, "
173 "const tmp<surfaceScalarField>&) : "
174 "interpolating volTypeField from cells to faces "
175 "without explicit correction"
179 const surfaceScalarField& lambdas = tlambdas();
180 const surfaceScalarField& ys = tys();
182 const Field<Type>& vfi = vf.internalField();
183 const scalarField& lambda = lambdas.internalField();
184 const scalarField& y = ys.internalField();
186 const fvMesh& mesh = vf.mesh();
187 const labelUList& P = mesh.owner();
188 const labelUList& N = mesh.neighbour();
190 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf
192 new GeometricField<Type, fvsPatchField, surfaceMesh>
196 "interpolate("+vf.name()+')',
204 GeometricField<Type, fvsPatchField, surfaceMesh>& sf = tsf();
206 Field<Type>& sfi = sf.internalField();
208 for (register label fi=0; fi<P.size(); fi++)
210 sfi[fi] = lambda[fi]*vfi[P[fi]] + y[fi]*vfi[N[fi]];
214 // Interpolate across coupled patches using given lambdas and ys
216 forAll(lambdas.boundaryField(), pi)
218 const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
219 const fvsPatchScalarField& pY = ys.boundaryField()[pi];
221 if (vf.boundaryField()[pi].coupled())
223 sf.boundaryField()[pi] =
224 pLambda*vf.boundaryField()[pi].patchInternalField()
225 + pY*vf.boundaryField()[pi].patchNeighbourField();
229 sf.boundaryField()[pi] = vf.boundaryField()[pi];
240 //- Return the face-interpolate of the given cell field
241 // with the given weigting factors
243 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
244 surfaceInterpolationScheme<Type>::interpolate
246 const GeometricField<Type, fvPatchField, volMesh>& vf,
247 const tmp<surfaceScalarField>& tlambdas
250 if (surfaceInterpolation::debug)
252 Info<< "surfaceInterpolationScheme<Type>::interpolate"
253 "(const GeometricField<Type, fvPatchField, volMesh>&, "
254 "const tmp<surfaceScalarField>&) : "
255 "interpolating volTypeField from cells to faces "
256 "without explicit correction"
260 const surfaceScalarField& lambdas = tlambdas();
262 const Field<Type>& vfi = vf.internalField();
263 const scalarField& lambda = lambdas.internalField();
265 const fvMesh& mesh = vf.mesh();
266 const labelUList& P = mesh.owner();
267 const labelUList& N = mesh.neighbour();
269 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf
271 new GeometricField<Type, fvsPatchField, surfaceMesh>
275 "interpolate("+vf.name()+')',
283 GeometricField<Type, fvsPatchField, surfaceMesh>& sf = tsf();
285 Field<Type>& sfi = sf.internalField();
287 for (register label fi=0; fi<P.size(); fi++)
289 sfi[fi] = lambda[fi]*(vfi[P[fi]] - vfi[N[fi]]) + vfi[N[fi]];
292 // Interpolate across coupled patches using given lambdas
294 forAll(lambdas.boundaryField(), pi)
296 const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
298 if (vf.boundaryField()[pi].coupled())
300 tsf().boundaryField()[pi] =
301 pLambda*vf.boundaryField()[pi].patchInternalField()
302 + (1.0 - pLambda)*vf.boundaryField()[pi].patchNeighbourField();
306 sf.boundaryField()[pi] = vf.boundaryField()[pi];
316 //- Return the face-interpolate of the given cell field
317 // with explicit correction
319 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
320 surfaceInterpolationScheme<Type>::interpolate
322 const GeometricField<Type, fvPatchField, volMesh>& vf
325 if (surfaceInterpolation::debug)
327 Info<< "surfaceInterpolationScheme<Type>::interpolate"
328 "(const GeometricField<Type, fvPatchField, volMesh>&) : "
329 << "interpolating volTypeField from cells to faces"
333 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf
334 = interpolate(vf, weights(vf));
338 tsf() += correction(vf);
345 //- Return the face-interpolate of the given cell field
346 // with explicit correction
348 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
349 surfaceInterpolationScheme<Type>::interpolate
351 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
354 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tinterpVf
355 = interpolate(tvf());
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
363 } // End namespace Foam
365 // ************************************************************************* //