surfaceInterpolationScheme.H
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 Class
27  Foam::surfaceInterpolationScheme
28 
29 Description
30  Abstract base class for surface interpolation schemes.
31 
32 SourceFiles
33  surfaceInterpolationScheme.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef surfaceInterpolationScheme_H
38 #define surfaceInterpolationScheme_H
39 
40 #include "tmp.H"
41 #include "volFieldsFwd.H"
42 #include "surfaceFieldsFwd.H"
43 #include "typeInfo.H"
44 #include "runTimeSelectionTables.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 class fvMesh;
52 
53 /*---------------------------------------------------------------------------*\
54  Class surfaceInterpolationScheme Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 template<class Type>
59 :
60  public refCount
61 {
62  // Private data
63 
64  //- Hold reference to mesh
65  const fvMesh& mesh_;
66 
67 
68  // Private Member Functions
69 
70  //- No copy construct
72 
73  //- No copy assignment
74  void operator=(const surfaceInterpolationScheme&) = delete;
75 
76 
77 public:
78 
79  //- Runtime type information
80  TypeName("surfaceInterpolationScheme");
81 
82 
83  // Declare run-time constructor selection tables
84 
86  (
87  tmp,
89  Mesh,
90  (
91  const fvMesh& mesh,
92  Istream& schemeData
93  ),
94  (mesh, schemeData)
95  );
96 
98  (
99  tmp,
101  MeshFlux,
102  (
103  const fvMesh& mesh,
104  const surfaceScalarField& faceFlux,
105  Istream& schemeData
106  ),
107  (mesh, faceFlux, schemeData)
108  );
109 
110 
111  // Constructors
112 
113  //- Construct from mesh
115  :
116  mesh_(mesh)
117  {}
118 
120  // Selectors
121 
122  //- Return new tmp interpolation scheme
124  (
125  const fvMesh& mesh,
126  Istream& schemeData
127  );
128 
129  //- Return new tmp interpolation scheme
131  (
132  const fvMesh& mesh,
133  const surfaceScalarField& faceFlux,
134  Istream& schemeData
135  );
136 
137 
138  //- Destructor
139  virtual ~surfaceInterpolationScheme() = default;
140 
141 
142  // Member Functions
143 
144  //- Return mesh reference
145  const fvMesh& mesh() const
146  {
147  return mesh_;
148  }
149 
150 
151  //- Return the face-interpolate of the given cell field
152  // with the given owner and neighbour weighting factors
155  (
159  );
160 
161  //- Return the face-interpolate of the given cell field
162  // with the given weighting factors dotted with given field Sf
163  template<class SFType>
164  static tmp
165  <
167  <
171  >
172  >
174  (
175  const SFType& Sf,
177  const tmp<surfaceScalarField>& tlambdas
178  );
179 
180  //- Return the face-interpolate of the given cell field
181  // with the given weighting factors
184  (
187  );
188 
189  //- Return the interpolation weighting factors for the given field
191  (
193  ) const = 0;
194 
195  //- Return true if this scheme uses an explicit correction
196  virtual bool corrected() const
197  {
198  return false;
199  }
200 
201  //- Return the explicit correction to the face-interpolate
202  // for the given field
205  {
206  return nullptr;
207  }
208 
209  //- Return the face-interpolate of the given cell field
210  // with explicit correction dotted with given field Sf
211  virtual
212  tmp
213  <
214  GeometricField
215  <
217  fvsPatchField,
218  surfaceMesh
219  >
220  >
222  (
223  const surfaceVectorField& Sf,
225  ) const;
226 
227  //- Return the face-interpolate of the given tmp cell field
228  // with explicit correction dotted with given field Sf
229  tmp
230  <
232  <
236  >
237  >
239  (
240  const surfaceVectorField& Sf,
242  ) const;
243 
244  //- Return the face-interpolate of the given cell field
245  // with explicit correction
248 
249  //- Return the face-interpolate of the given tmp cell field
250  // with explicit correction
253  (
255  ) const;
256 };
257 
258 
259 template<>
260 tmp
261 <
263  <
267  >
268 >
270 (
271  const surfaceVectorField& Sf,
273 ) const;
274 
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace Foam
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 // Add the patch constructor functions to the hash tables
283 
284 #define makeSurfaceInterpolationTypeScheme(SS, Type) \
285  \
286 defineNamedTemplateTypeNameAndDebug(SS<Type>, 0); \
287  \
288 surfaceInterpolationScheme<Type>::addMeshConstructorToTable<SS<Type>> \
289  add##SS##Type##MeshConstructorToTable_; \
290  \
291 surfaceInterpolationScheme<Type>::addMeshFluxConstructorToTable<SS<Type>> \
292  add##SS##Type##MeshFluxConstructorToTable_;
293 
294 #define makeSurfaceInterpolationScheme(SS) \
295  \
296 makeSurfaceInterpolationTypeScheme(SS, scalar) \
297 makeSurfaceInterpolationTypeScheme(SS, vector) \
298 makeSurfaceInterpolationTypeScheme(SS, sphericalTensor) \
299 makeSurfaceInterpolationTypeScheme(SS, symmTensor) \
300 makeSurfaceInterpolationTypeScheme(SS, tensor)
301 
302 
303 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304 
305 #ifdef NoRepository
307 #endif
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 #endif
312 
313 // ************************************************************************* //
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
declareRunTimeSelectionTable(tmp, surfaceInterpolationScheme, Mesh,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:43
Reference counter for various OpenFOAM components.
Definition: refCount.H:44
Forwards and collection of common volume field types.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
virtual ~surfaceInterpolationScheme()=default
Destructor.
static tmp< GeometricField< typename innerProduct< typename SFType::value_type, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const SFType &Sf, const GeometricField< Type, fvPatchField, volMesh > &vf, const tmp< surfaceScalarField > &tlambdas)
Return the face-interpolate of the given cell field.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors for the given field.
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
TypeName("surfaceInterpolationScheme")
Runtime type information.
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) - 2 >::type type
Definition: products.H:155
Abstract base class for surface interpolation schemes.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
Namespace for OpenFOAM.