surfaceInterpolate.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
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 InNamespace
25  Foam::fvc
26 
27 Description
28  Surface Interpolation
29 
30 SourceFiles
31  surfaceInterpolate.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef surfaceInterpolate_H
36 #define surfaceInterpolate_H
37 
38 #include "tmp.H"
39 #include "volFieldsFwd.H"
40 #include "surfaceFieldsFwd.H"
42 #include "one.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Namespace fvc functions Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 namespace fvc
54 {
55  //- Return weighting factors for scheme given from Istream
56  template<class Type>
57  static tmp<surfaceInterpolationScheme<Type>> scheme
58  (
59  const surfaceScalarField& faceFlux,
60  Istream& schemeData
61  );
62 
63  //- Return weighting factors for scheme given by name in dictionary
64  template<class Type>
65  static tmp<surfaceInterpolationScheme<Type>> scheme
66  (
67  const surfaceScalarField& faceFlux,
68  const word& name
69  );
70 
71 
72  //- Return weighting factors for scheme given from Istream
73  template<class Type>
74  static tmp<surfaceInterpolationScheme<Type>> scheme
75  (
76  const fvMesh& mesh,
77  Istream& schemeData
78  );
79 
80  //- Return weighting factors for scheme given by name in dictionary
81  template<class Type>
82  static tmp<surfaceInterpolationScheme<Type>> scheme
83  (
84  const fvMesh& mesh,
85  const word& name
86  );
87 
88 
89  //- Interpolate field onto faces using scheme given by Istream
90  template<class Type>
91  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
92  (
93  const GeometricField<Type, fvPatchField, volMesh>& tvf,
94  const surfaceScalarField& faceFlux,
95  Istream& schemeData
96  );
97 
98  //- Interpolate field onto faces using scheme given by name in fvSchemes
99  template<class Type>
100  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
101  (
102  const GeometricField<Type, fvPatchField, volMesh>& tvf,
103  const surfaceScalarField& faceFlux,
104  const word& name
105  );
106 
107  //- Interpolate field onto faces using scheme given by name in fvSchemes
108  template<class Type>
109  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
110  (
111  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf,
112  const surfaceScalarField& faceFlux,
113  const word& name
114  );
115 
116  //- Interpolate field onto faces using scheme given by name in fvSchemes
117  template<class Type>
118  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
119  (
120  const GeometricField<Type, fvPatchField, volMesh>& tvf,
121  const tmp<surfaceScalarField>& faceFlux,
122  const word& name
123  );
124 
125  //- Interpolate field onto faces using scheme given by name in fvSchemes
126  template<class Type>
127  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
128  (
129  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf,
130  const tmp<surfaceScalarField>& faceFlux,
131  const word& name
132  );
133 
134 
135  //- Interpolate field onto faces using scheme given by Istream
136  template<class Type>
137  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
138  (
139  const GeometricField<Type, fvPatchField, volMesh>& tvf,
140  Istream& schemeData
141  );
142 
143  //- Interpolate field onto faces using scheme given by name in fvSchemes
144  template<class Type>
145  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
146  (
147  const GeometricField<Type, fvPatchField, volMesh>& tvf,
148  const word& name
149  );
150 
151  //- Interpolate field onto faces using scheme given by name in fvSchemes
152  template<class Type>
153  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
154  (
155  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf,
156  const word& name
157  );
158 
159 
160  //- Interpolate field onto faces using 'interpolate(<name>)'
161  template<class Type>
162  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
163  (
164  const GeometricField<Type, fvPatchField, volMesh>& tvf
165  );
166 
167  //- Interpolate tmp field onto faces using 'interpolate(<name>)'
168  template<class Type>
169  static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> interpolate
170  (
171  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf
172  );
173 
174 
175  //- Interpolate boundary field onto faces (simply a type conversion)
176  template<class Type>
177  static tmp<FieldField<fvsPatchField, Type>> interpolate
178  (
179  const FieldField<fvPatchField, Type>& fvpff
180  );
181 
182  //- Interpolate boundary field onto faces (simply a type conversion)
183  template<class Type>
184  static tmp<FieldField<fvsPatchField, Type>> interpolate
185  (
186  const tmp<FieldField<fvPatchField, Type>>& tfvpff
187  );
188 
189  //- Interpolate 'one' returning 'one'
190  inline one interpolate(const one&)
191  {
192  return one();
193  }
194 
195 
196  //- Interpolate field onto faces
197  // and 'dot' with given surfaceVectorField Sf
198  template<class Type>
199  static
200  tmp
201  <
203  <
207  >
209  (
210  const surfaceVectorField& Sf,
212  );
213 
214  //- Interpolate tmp field onto faces
215  // and 'dot' with given surfaceVectorField Sf
216  template<class Type>
217  static
218  tmp
219  <
221  <
222  typename innerProduct<vector, Type>::type,
223  fvsPatchField,
225  >
227  (
228  const surfaceVectorField& Sf,
230  );
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 } // End namespace Foam
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #ifdef NoRepository
241  #include "surfaceInterpolate.C"
242 #endif
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 #endif
247 
248 // ************************************************************************* //
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:47
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
static tmp< surfaceInterpolationScheme< Type > > scheme(const surfaceScalarField &faceFlux, Istream &schemeData)
Return weighting factors for scheme given from Istream.
Generic GeometricField class.
dynamicFvMesh & mesh
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:50