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-2022 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 "MeshObject.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 :
56  public MeshObject<fvMesh, UpdateableMeshObject, volPointInterpolation>
57 {
58  // Private Data
59 
60  //- Interpolation scheme weighting factor array.
61  scalarListList pointWeights_;
62 
63  //- Boundary addressing
64  autoPtr<primitivePatch> boundaryPtr_;
65 
66  //- Per boundary point the weights per pointFaces.
67  List<scalarListList> boundaryPointWeights_;
68 
69  //- Per boundary point the weights per pointFaces.
70  List<scalarListList> boundaryPointNbrWeights_;
71 
72 
73  // Private Member Functions
74 
75  //- Construct all point weighting factors
76  void makeWeights();
77 
78 
79 public:
80 
81  // Declare name of the class and its debug switch
82  ClassName("volPointInterpolation");
83 
84 
85  // Constructors
86 
87  //- Constructor given fvMesh and pointMesh.
88  explicit volPointInterpolation(const fvMesh&);
89 
90  //- Disallow default bitwise copy construction
92 
93 
94  //- Destructor
96 
97 
98  // Member Functions
99 
100  //- Correct weighting factors for moving mesh.
101  virtual bool movePoints();
102 
103  //- Update mesh topology using the morph engine
104  virtual void topoChange(const polyTopoChangeMap&);
105 
106  //- Update from another mesh using the given map
107  virtual void mapMesh(const polyMeshMap&);
108 
109  //- Redistribute or update using the given distribution map
110  virtual void distribute(const polyDistributionMap&);
111 
112 
113  // Interpolation functions
114 
115  //- Interpolate volField using inverse distance weighting
116  // returning pointField
117  template<class Type>
119  (
121  ) const;
122 
123  //- Interpolate tmp<volField> using inverse distance weighting
124  // returning pointField
125  template<class Type>
127  (
129  ) const;
130 
131 
132  // Low level
133 
134  //- Interpolate from volField to pointField
135  // using inverse distance weighting
136  template<class Type>
138  (
141  ) const;
142 
143  //- Interpolate from volField to pointField
144  // using inverse distance weighting
145  template<class Type>
146  void interpolate
147  (
150  ) const;
151 
152  //- Interpolate volField using inverse distance weighting
153  // returning pointField with name. Optionally caches.
154  template<class Type>
156  (
158  const word& name,
159  const bool cache
160  ) const;
161 
162 
163  // Interpolation for displacement (applies 2D correction)
164 
165  //- Interpolate from volField to pointField
166  // using inverse distance weighting
168  (
169  const volVectorField&,
171  ) const;
172 
173 
174  // Member Operators
175 
176  //- Disallow default bitwise assignment
177  void operator=(const volPointInterpolation&) = delete;
178 };
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 } // End namespace Foam
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 #ifdef NoRepository
189 #endif
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 #endif
194 
195 // ************************************************************************* //
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
const word & name() const
Return name.
Definition: IOobject.H:315
Generic GeometricField class.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
ClassName("volPointInterpolation")
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:87
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
void interpolateUnconstrained(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate from volField to pointField.
virtual void topoChange(const polyTopoChangeMap &)
Update mesh topology using the morph engine.
Interpolate from cell centres to points (vertices) using inverse distance weighting.
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
volPointInterpolation(const fvMesh &)
Constructor given fvMesh and pointMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
virtual bool movePoints()
Correct weighting factors for moving mesh.
Namespace for OpenFOAM.
void operator=(const volPointInterpolation &)=delete
Disallow default bitwise assignment.