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 "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  UpdateableMeshObject,
60  volPointInterpolation
61  >
62 {
63  // Private Data
64 
65  //- Interpolation scheme weighting factor array.
66  scalarListList pointWeights_;
67 
68  //- Boundary addressing
69  autoPtr<primitivePatch> boundaryPtr_;
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  // Private Member Functions
79 
80  //- Construct all point weighting factors
81  void makeWeights();
82 
83 
84 protected:
85 
86  friend class DemandDrivenMeshObject
87  <
88  fvMesh,
91  >;
92 
93  // Protected Constructors
94 
95  //- Constructor given fvMesh
96  explicit volPointInterpolation(const fvMesh&);
97 
98 
99 public:
100 
101  // Declare name of the class and its debug switch
102  ClassName("volPointInterpolation");
103 
104 
105  // Constructors
106 
107  //- Disallow default bitwise copy construction
109 
110 
111  //- Destructor
113 
114 
115  // Member Functions
116 
117  //- Correct weighting factors for moving mesh.
118  virtual bool movePoints();
119 
120  //- Update mesh topology using the morph engine
121  virtual void topoChange(const polyTopoChangeMap&);
122 
123  //- Update from another mesh using the given map
124  virtual void mapMesh(const polyMeshMap&);
125 
126  //- Redistribute or update using the given distribution map
127  virtual void distribute(const polyDistributionMap&);
128 
129 
130  // Interpolation functions
131 
132  //- Interpolate volField using inverse distance weighting
133  // returning pointField
134  template<class Type>
136  (
137  const VolField<Type>&
138  ) const;
139 
140  //- Interpolate tmp<volField> using inverse distance weighting
141  // returning pointField
142  template<class Type>
144  (
145  const tmp<VolField<Type>>&
146  ) const;
147 
148 
149  // Low level
150 
151  //- Interpolate from volField to pointField
152  // using inverse distance weighting
153  template<class Type>
155  (
156  const VolField<Type>&,
158  ) const;
159 
160  //- Interpolate from volField to pointField
161  // using inverse distance weighting
162  template<class Type>
163  void interpolate
164  (
165  const VolField<Type>&,
167  ) const;
168 
169  //- Interpolate volField using inverse distance weighting
170  // returning pointField with name. Optionally caches.
171  template<class Type>
173  (
174  const VolField<Type>&,
175  const word& name,
176  const bool cache
177  ) const;
178 
179 
180  // Interpolation for displacement (applies 2D correction)
181 
182  //- Interpolate from volField to pointField
183  // using inverse distance weighting
185  (
186  const volVectorField&,
188  ) const;
189 
190 
191  // Member Operators
192 
193  //- Disallow default bitwise assignment
194  void operator=(const volPointInterpolation&) = delete;
195 };
196 
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 } // End namespace Foam
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 #ifdef NoRepository
206 #endif
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 #endif
211 
212 // ************************************************************************* //
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
Interpolate from cell centres to points (vertices) using inverse distance weighting.
virtual bool movePoints()
Correct weighting factors for moving mesh.
void operator=(const volPointInterpolation &)=delete
Disallow default bitwise assignment.
virtual void topoChange(const polyTopoChangeMap &)
Update mesh topology using the morph engine.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
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.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
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.