Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
finiteVolume
interpolation
surfaceInterpolation
surfaceInterpolation
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 | Copyright (C) 1991-2010 OpenCFD Ltd.
6
\\/ M anipulation |
7
-------------------------------------------------------------------------------
8
License
9
This file is part of OpenFOAM.
10
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.
15
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.
20
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/>.
23
24
\*---------------------------------------------------------------------------*/
25
26
#include <
finiteVolume/surfaceInterpolate.H
>
27
28
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29
30
namespace
Foam
31
{
32
33
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35
namespace
fvc
36
{
37
38
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40
// Return weighting factors for scheme given by name in dictionary
41
template
<
class
Type>
42
tmp<surfaceInterpolationScheme<Type> >
scheme
43
(
44
const
surfaceScalarField
& faceFlux,
45
Istream
& streamData
46
)
47
{
48
return
surfaceInterpolationScheme<Type>::New
49
(
50
faceFlux.
mesh
(),
51
faceFlux,
52
streamData
53
);
54
}
55
56
57
// Return weighting factors for scheme given by name in dictionary
58
template
<
class
Type>
59
tmp<surfaceInterpolationScheme<Type>
>
scheme
60
(
61
const
surfaceScalarField
& faceFlux,
62
const
word
&
name
63
)
64
{
65
return
surfaceInterpolationScheme<Type>::New
66
(
67
faceFlux.
mesh
(),
68
faceFlux,
69
faceFlux.
mesh
().
interpolationScheme
(name)
70
);
71
}
72
73
74
// Return weighting factors for scheme given by name in dictionary
75
template
<
class
Type>
76
tmp<surfaceInterpolationScheme<Type>
>
scheme
77
(
78
const
fvMesh
&
mesh
,
79
Istream
& streamData
80
)
81
{
82
return
surfaceInterpolationScheme<Type>::New
83
(
84
mesh,
85
streamData
86
);
87
}
88
89
90
// Return weighting factors for scheme given by name in dictionary
91
template
<
class
Type>
92
tmp<surfaceInterpolationScheme<Type>
>
scheme
93
(
94
const
fvMesh
& mesh,
95
const
word
& name
96
)
97
{
98
return
surfaceInterpolationScheme<Type>::New
99
(
100
mesh,
101
mesh.
interpolationScheme
(name)
102
);
103
}
104
105
106
// Interpolate field onto faces using scheme given by name in dictionary
107
template
<
class
Type>
108
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
109
interpolate
110
(
111
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
112
const
surfaceScalarField
& faceFlux,
113
Istream
& schemeData
114
)
115
{
116
if
(
surfaceInterpolation::debug
)
117
{
118
Info
<<
"interpolate"
119
<<
"(const GeometricField<Type, fvPatchField, volMesh>&, "
120
<<
"const surfaceScalarField&, Istream&) : "
121
<<
"interpolating GeometricField<Type, fvPatchField, volMesh> "
122
<<
endl
;
123
}
124
125
return
scheme<Type>(faceFlux, schemeData)().interpolate(vf);
126
}
127
128
129
// Interpolate field onto faces using scheme given by name in dictionary
130
template
<
class
Type>
131
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
132
interpolate
133
(
134
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
135
const
surfaceScalarField
& faceFlux,
136
const
word
& name
137
)
138
{
139
if
(
surfaceInterpolation::debug
)
140
{
141
Info
<<
"interpolate"
142
<<
"(const GeometricField<Type, fvPatchField, volMesh>&, "
143
<<
"const surfaceScalarField&, const word&) : "
144
<<
"interpolating GeometricField<Type, fvPatchField, volMesh> "
145
<<
"using "
<< name
146
<<
endl
;
147
}
148
149
return
scheme<Type>(faceFlux,
name
)().interpolate(vf);
150
}
151
152
// Interpolate field onto faces using scheme given by name in dictionary
153
template
<
class
Type>
154
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
155
interpolate
156
(
157
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
158
const
surfaceScalarField
& faceFlux,
159
const
word
& name
160
)
161
{
162
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
> tsf =
163
interpolate
(tvf(), faceFlux, name);
164
165
tvf.clear();
166
167
return
tsf;
168
}
169
170
// Interpolate field onto faces using scheme given by name in dictionary
171
template
<
class
Type>
172
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
173
interpolate
174
(
175
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
176
const
tmp<surfaceScalarField>
& tFaceFlux,
177
const
word
& name
178
)
179
{
180
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
> tsf =
181
interpolate
(vf, tFaceFlux(), name);
182
183
tFaceFlux.
clear
();
184
185
return
tsf;
186
}
187
188
// Interpolate field onto faces using scheme given by name in dictionary
189
template
<
class
Type>
190
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
191
interpolate
192
(
193
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
194
const
tmp<surfaceScalarField>
& tFaceFlux,
195
const
word
& name
196
)
197
{
198
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
> tsf =
199
interpolate
(tvf(), tFaceFlux(), name);
200
201
tvf.clear();
202
tFaceFlux.
clear
();
203
204
return
tsf;
205
}
206
207
208
// Interpolate field onto faces using scheme given by name in dictionary
209
template
<
class
Type>
210
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
211
interpolate
212
(
213
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
214
Istream
& schemeData
215
)
216
{
217
if
(
surfaceInterpolation::debug
)
218
{
219
Info
<<
"interpolate"
220
<<
"(const GeometricField<Type, fvPatchField, volMesh>&, "
221
<<
"Istream&) : "
222
<<
"interpolating GeometricField<Type, fvPatchField, volMesh> "
223
<<
endl
;
224
}
225
226
return
scheme<Type>(vf.
mesh
(), schemeData)().interpolate(vf);
227
}
228
229
// Interpolate field onto faces using scheme given by name in dictionary
230
template
<
class
Type>
231
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
232
interpolate
233
(
234
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
235
const
word
& name
236
)
237
{
238
if
(
surfaceInterpolation::debug
)
239
{
240
Info
<<
"interpolate"
241
<<
"(const GeometricField<Type, fvPatchField, volMesh>&, "
242
<<
"const word&) : "
243
<<
"interpolating GeometricField<Type, fvPatchField, volMesh> "
244
<<
"using "
<< name
245
<<
endl
;
246
}
247
248
return
scheme<Type>(vf.
mesh
(),
name
)().interpolate(vf);
249
}
250
251
// Interpolate field onto faces using scheme given by name in dictionary
252
template
<
class
Type>
253
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
254
interpolate
255
(
256
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
257
const
word
& name
258
)
259
{
260
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
> tsf =
261
interpolate
(tvf(), name);
262
263
tvf.clear();
264
265
return
tsf;
266
}
267
268
269
// Interpolate field onto faces using central differencing
270
template
<
class
Type>
271
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
272
interpolate
273
(
274
const
GeometricField<Type, fvPatchField, volMesh>
& vf
275
)
276
{
277
if
(
surfaceInterpolation::debug
)
278
{
279
Info
<<
"interpolate"
280
<<
"(const GeometricField<Type, fvPatchField, volMesh>&) : "
281
<<
"interpolating GeometricField<Type, fvPatchField, volMesh> "
282
<<
"using run-time selected scheme"
283
<<
endl
;
284
}
285
286
return
interpolate
(vf,
"interpolate("
+ vf.
name
() +
')'
);
287
}
288
289
290
// Interpolate field onto faces using central differencing
291
template
<
class
Type>
292
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
>
293
interpolate
294
(
295
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf
296
)
297
{
298
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
> tsf =
299
interpolate
(tvf());
300
tvf.clear();
301
return
tsf;
302
}
303
304
305
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306
307
}
// End namespace fvc
308
309
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310
311
}
// End namespace Foam
312
313
// ************************ vim: set sw=4 sts=4 et: ************************ //