surfaceInterpolate.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "surfaceInterpolate.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 template<class Type>
35 (
36  const surfaceScalarField& faceFlux,
37  Istream& streamData
38 )
39 {
41  (
42  faceFlux.mesh(),
43  faceFlux,
44  streamData
45  );
46 }
47 
48 
49 template<class Type>
51 (
52  const surfaceScalarField& faceFlux,
53  const word& name
54 )
55 {
57  (
58  faceFlux.mesh(),
59  faceFlux,
60  faceFlux.mesh().interpolationScheme(name)
61  );
62 }
63 
64 
65 template<class Type>
67 (
68  const fvMesh& mesh,
69  Istream& streamData
70 )
71 {
73  (
74  mesh,
75  streamData
76  );
77 }
78 
79 
80 template<class Type>
82 (
83  const fvMesh& mesh,
84  const word& name
85 )
86 {
88  (
89  mesh,
90  mesh.interpolationScheme(name)
91  );
92 }
93 
94 
95 template<class Type>
98 (
99  const GeometricField<Type, fvPatchField, volMesh>& vf,
100  const surfaceScalarField& faceFlux,
101  Istream& schemeData
102 )
103 {
105  {
107  << "interpolating GeometricField<Type, fvPatchField, volMesh> "
108  << vf.name() << endl;
109  }
110 
111  return scheme<Type>(faceFlux, schemeData)().interpolate(vf);
112 }
113 
114 
115 template<class Type>
118 (
119  const GeometricField<Type, fvPatchField, volMesh>& vf,
120  const surfaceScalarField& faceFlux,
121  const word& name
122 )
123 {
125  {
127  << "interpolating GeometricField<Type, fvPatchField, volMesh> "
128  << vf.name() << " using " << name << endl;
129  }
130 
131  return scheme<Type>(faceFlux, name)().interpolate(vf);
132 }
133 
134 template<class Type>
137 (
138  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf,
139  const surfaceScalarField& faceFlux,
140  const word& name
141 )
142 {
143  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf =
144  interpolate(tvf(), faceFlux, name);
145 
146  tvf.clear();
147 
148  return tsf;
149 }
150 
151 template<class Type>
154 (
155  const GeometricField<Type, fvPatchField, volMesh>& vf,
156  const tmp<surfaceScalarField>& tFaceFlux,
157  const word& name
158 )
159 {
160  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf =
161  interpolate(vf, tFaceFlux(), name);
162 
163  tFaceFlux.clear();
164 
165  return tsf;
166 }
167 
168 template<class Type>
171 (
172  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf,
173  const tmp<surfaceScalarField>& tFaceFlux,
174  const word& name
175 )
176 {
177  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf =
178  interpolate(tvf(), tFaceFlux(), name);
179 
180  tvf.clear();
181  tFaceFlux.clear();
182 
183  return tsf;
184 }
185 
186 
187 template<class Type>
190 (
191  const GeometricField<Type, fvPatchField, volMesh>& vf,
192  Istream& schemeData
193 )
194 {
196  {
198  << "interpolating GeometricField<Type, fvPatchField, volMesh> "
199  << vf.name() << endl;
200  }
201 
202  return scheme<Type>(vf.mesh(), schemeData)().interpolate(vf);
203 }
204 
205 template<class Type>
208 (
209  const GeometricField<Type, fvPatchField, volMesh>& vf,
210  const word& name
211 )
212 {
214  {
216  << "interpolating GeometricField<Type, fvPatchField, volMesh> "
217  << vf.name() << " using " << name
218  << endl;
219  }
220 
221  return scheme<Type>(vf.mesh(), name)().interpolate(vf);
222 }
223 
224 template<class Type>
227 (
228  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf,
229  const word& name
230 )
231 {
232  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf =
233  interpolate(tvf(), name);
234 
235  tvf.clear();
236 
237  return tsf;
238 }
239 
240 
241 template<class Type>
244 (
245  const GeometricField<Type, fvPatchField, volMesh>& vf
246 )
247 {
249  {
251  << "interpolating GeometricField<Type, fvPatchField, volMesh> "
252  << vf.name() << " using run-time selected scheme"
253  << endl;
254  }
255 
256  return interpolate(vf, "interpolate(" + vf.name() + ')');
257 }
258 
259 
260 template<class Type>
263 (
264  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf
265 )
266 {
267  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf =
268  interpolate(tvf());
269  tvf.clear();
270  return tsf;
271 }
272 
273 
274 template<class Type>
277 (
278  const FieldField<fvPatchField, Type>& fvpff
279 )
280 {
281  auto tresult = tmp<FieldField<fvsPatchField, Type>>::New(fvpff.size());
282  auto& result = tresult.ref();
283 
284  forAll(result, patchi)
285  {
286  result.set
287  (
288  patchi,
289  fvsPatchField<Type>::NewCalculatedType(fvpff[patchi].patch()).ptr()
290  );
291  result[patchi] = fvpff[patchi];
292  }
293 
294  return tresult;
295 }
296 
297 
298 template<class Type>
301 (
302  const tmp<FieldField<fvPatchField, Type>>& tfvpff
303 )
304 {
305  tmp<FieldField<fvsPatchField, Type>> tfvspff = interpolate(tfvpff());
306  tfvpff.clear();
307  return tfvspff;
308 }
309 
310 
311 template<class Type>
312 Foam::tmp
313 <
315  <
319  >
320 >
322 (
323  const surfaceVectorField& Sf,
324  const GeometricField<Type, fvPatchField, volMesh>& vf
325 )
326 {
328  {
330  << "interpolating GeometricField<Type, fvPatchField, volMesh> "
331  << vf.name() << " using run-time selected scheme"
332  << endl;
333  }
334 
335  return scheme<Type>
336  (
337  vf.mesh(),
338  "dotInterpolate(" + Sf.name() + ',' + vf.name() + ')'
339  )().dotInterpolate(Sf, vf);
340 }
341 
342 
343 template<class Type>
344 Foam::tmp
345 <
347  <
351  >
352 >
354 (
355  const surfaceVectorField& Sf,
356  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf
357 )
358 {
359  tmp
360  <
361  GeometricField
362  <
364  fvsPatchField,
365  surfaceMesh
366  >
367  > tsf = dotInterpolate(Sf, tvf());
368  tvf.clear();
369  return tsf;
370 }
371 
372 
373 // ************************************************************************* //
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:43
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
static tmp< surfaceInterpolationScheme< Type > > scheme(const surfaceScalarField &faceFlux, Istream &schemeData)
Return weighting factors for scheme given from Istream.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Generic GeometricField class.
Surface Interpolation.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
int debug
Static debugging option.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
const std::string patch
OpenFOAM patch number as a std::string.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) - 2 >::type type
Definition: products.H:155
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
#define InfoInFunction
Report an information message using Foam::Info.