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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::surfaceInterpolationScheme
28
29Description
30 Abstract base class for surface interpolation schemes.
31
32SourceFiles
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"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51class fvMesh;
52
53/*---------------------------------------------------------------------------*\
54 Class surfaceInterpolationScheme Declaration
55\*---------------------------------------------------------------------------*/
56
57template<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
77public:
78
79 //- Runtime type information
80 TypeName("surfaceInterpolationScheme");
81
82
83 // Declare run-time constructor selection tables
86 (
87 tmp,
89 Mesh,
90 (
91 const fvMesh& mesh,
92 Istream& schemeData
93 ),
94 (mesh, schemeData)
95 );
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
119
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 >
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 {
207 (
208 nullptr
209 );
210 }
211
212 //- Return the face-interpolate of the given cell field
213 // with explicit correction dotted with given field Sf
214 virtual
215 tmp
216 <
218 <
222 >
223 >
225 (
226 const surfaceVectorField& Sf,
228 ) const;
229
230 //- Return the face-interpolate of the given tmp cell field
231 // with explicit correction dotted with given field Sf
232 tmp
233 <
235 <
239 >
240 >
242 (
243 const surfaceVectorField& Sf,
245 ) const;
246
247 //- Return the face-interpolate of the given cell field
248 // with explicit correction
251
252 //- Return the face-interpolate of the given tmp cell field
253 // with explicit correction
256 (
258 ) const;
259};
260
261
262template<>
263tmp
264<
265 GeometricField
266 <
268 fvsPatchField,
269 surfaceMesh
270 >
273(
274 const surfaceVectorField& Sf,
276) const;
277
278
279// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280
281} // End namespace Foam
282
283// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284
285// Add the patch constructor functions to the hash tables
287#define makeSurfaceInterpolationTypeScheme(SS, Type) \
288 \
289defineNamedTemplateTypeNameAndDebug(SS<Type>, 0); \
290 \
291surfaceInterpolationScheme<Type>::addMeshConstructorToTable<SS<Type>> \
292 add##SS##Type##MeshConstructorToTable_; \
293 \
294surfaceInterpolationScheme<Type>::addMeshFluxConstructorToTable<SS<Type>> \
295 add##SS##Type##MeshFluxConstructorToTable_;
297#define makeSurfaceInterpolationScheme(SS) \
298 \
299makeSurfaceInterpolationTypeScheme(SS, scalar) \
300makeSurfaceInterpolationTypeScheme(SS, vector) \
301makeSurfaceInterpolationTypeScheme(SS, sphericalTensor) \
302makeSurfaceInterpolationTypeScheme(SS, symmTensor) \
303makeSurfaceInterpolationTypeScheme(SS, tensor)
304
305
306// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307
308#ifdef NoRepository
310#endif
311
312// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313
314#endif
315
316// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) -2 >::type type
Definition: products.H:149
Reference counter for various OpenFOAM components.
Definition: refCount.H:51
Abstract base class for surface interpolation schemes.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
surfaceInterpolationScheme(const fvMesh &mesh)
Construct from mesh.
virtual tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors for the given field.
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.
TypeName("surfaceInterpolationScheme")
Runtime type information.
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.
declareRunTimeSelectionTable(tmp, surfaceInterpolationScheme, MeshFlux,(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData),(mesh, faceFlux, schemeData))
const fvMesh & mesh() const
Return mesh reference.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
declareRunTimeSelectionTable(tmp, surfaceInterpolationScheme, Mesh,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
virtual ~surfaceInterpolationScheme()=default
Destructor.
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:52
A class for managing temporary objects.
Definition: tmp.H:65
type
Volume classification types.
Definition: volumeType.H:66
Namespace for OpenFOAM.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73