volPointInterpolation.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-2023 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 Class
25  Foam::volPointInterpolation
26 
27 Description
28  Interpolate from cell centres to points (vertices) using inverse distance
29  weighting
30 
31 SourceFiles
32  volPointInterpolation.C
33  volPointInterpolationTemplates.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef volPointInterpolation_H
38 #define volPointInterpolation_H
39 
40 #include "DemandDrivenMeshObject.H"
41 #include "scalarList.H"
42 #include "volFields.H"
43 #include "pointFields.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class volPointInterpolation Declaration
52 \*---------------------------------------------------------------------------*/
53 
55 :
57  <
58  fvMesh,
59  DeletableMeshObject,
60  volPointInterpolation
61  >
62 {
63  // Private Data
64 
65  //- Interpolation scheme weighting factor array.
66  scalarListList pointWeights_;
67 
68  //- Boundary addressing
69  primitivePatch boundary_;
70 
71  //- Per boundary point the weights per pointFaces.
72  List<scalarListList> boundaryPointWeights_;
73 
74  //- Per boundary point the weights per pointFaces.
75  List<scalarListList> boundaryPointNbrWeights_;
76 
77 
78 protected:
79 
80  friend class DemandDrivenMeshObject
81  <
82  fvMesh,
85  >;
86 
87  // Protected Constructors
88 
89  //- Constructor given fvMesh
90  explicit volPointInterpolation(const fvMesh&);
91 
92 
93 public:
94 
95  // Declare name of the class and its debug switch
96  ClassName("volPointInterpolation");
97 
98 
99  // Constructors
100 
101  //- Disallow default bitwise copy construction
103 
104 
105  //- Destructor
107 
108 
109  // Member Functions
110 
111  //- Interpolate volField using inverse distance weighting
112  // returning pointField
113  template<class Type>
115  (
116  const VolField<Type>&
117  ) const;
118 
119  //- Interpolate tmp<volField> using inverse distance weighting
120  // returning pointField
121  template<class Type>
123  (
124  const tmp<VolField<Type>>&
125  ) const;
126 
127 
128  // Low level
129 
130  //- Interpolate from volField to pointField
131  // using inverse distance weighting
132  template<class Type>
134  (
135  const VolField<Type>&,
137  ) const;
138 
139  //- Interpolate from volField to pointField
140  // using inverse distance weighting
141  template<class Type>
142  void interpolate
143  (
144  const VolField<Type>&,
146  ) const;
147 
148  //- Interpolate volField using inverse distance weighting
149  // returning pointField with name. Optionally caches.
150  template<class Type>
152  (
153  const VolField<Type>&,
154  const word& name,
155  const bool cache
156  ) const;
157 
158 
159  // Interpolation for displacement (applies 2D correction)
160 
161  //- Interpolate from volField to pointField
162  // using inverse distance weighting
164  (
165  const volVectorField&,
167  ) const;
168 
169 
170  // Member Operators
171 
172  //- Disallow default bitwise assignment
173  void operator=(const volPointInterpolation&) = delete;
174 };
175 
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace Foam
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 #ifdef NoRepository
185 #endif
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
MeshObject types:
Definition: MeshObjects.H:67
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
A class for managing temporary objects.
Definition: tmp.H:55
Interpolate from cell centres to points (vertices) using inverse distance weighting.
void operator=(const volPointInterpolation &)=delete
Disallow default bitwise assignment.
volPointInterpolation(const fvMesh &)
Constructor given fvMesh.
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
void interpolateUnconstrained(const VolField< Type > &, PointField< Type > &) const
Interpolate from volField to pointField.
ClassName("volPointInterpolation")
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.