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-2018 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  volPointInterpolate.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 class fvMesh;
51 class pointMesh;
52 
53 /*---------------------------------------------------------------------------*\
54  Class volPointInterpolation Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 :
59  public MeshObject<fvMesh, UpdateableMeshObject, volPointInterpolation>
60 {
61  // Private data
62 
63  //- Interpolation scheme weighting factor array.
64  scalarListList pointWeights_;
65 
66 
67  // Boundary handling
68 
69  //- Boundary addressing
70  autoPtr<primitivePatch> boundaryPtr_;
71 
72  //- Per boundary face whether is on non-coupled patch
73  boolList boundaryIsPatchFace_;
74 
75  //- Per mesh(!) point whether is on non-coupled patch (on any
76  // processor)
77  boolList isPatchPoint_;
78 
79  //- Per boundary point the weights per pointFaces.
80  scalarListList boundaryPointWeights_;
81 
82 
83  // Private Member Functions
84 
85  //- Construct addressing over all boundary faces
86  void calcBoundaryAddressing();
87 
88  //- Make weights for internal and coupled-only boundarypoints
89  void makeInternalWeights(scalarField& sumWeights);
90 
91  //- Make weights for points on uncoupled patches
92  void makeBoundaryWeights(scalarField& sumWeights);
93 
94  //- Construct all point weighting factors
95  void makeWeights();
96 
97  //- Helper: push master point data to collocated points
98  template<class Type>
99  void pushUntransformedData(List<Type>&) const;
100 
101  //- Get boundary field in same order as boundary faces. Field is
102  // zero on all coupled and empty patches
103  template<class Type>
104  tmp<Field<Type>> flatBoundaryField
105  (
107  ) const;
108 
109  //- Add separated contributions
110  template<class Type>
111  void addSeparated
112  (
114  ) const;
115 
116  //- Disallow default bitwise copy construct
118 
119  //- Disallow default bitwise assignment
120  void operator=(const volPointInterpolation&);
121 
122 
123 public:
124 
125  // Declare name of the class and its debug switch
126  ClassName("volPointInterpolation");
127 
128 
129  // Constructors
130 
131  //- Constructor given fvMesh and pointMesh.
132  explicit volPointInterpolation(const fvMesh&);
133 
134 
135  //- Destructor
137 
138 
139  // Member functions
140 
141  // Edit
142 
143  //- Update mesh topology using the morph engine
144  void updateMesh(const mapPolyMesh&);
145 
146  //- Correct weighting factors for moving mesh.
147  bool movePoints();
148 
149 
150  // Interpolation functions
151 
152  //- Interpolate volField using inverse distance weighting
153  // returning pointField
154  template<class Type>
156  (
158  ) const;
159 
160  //- Interpolate tmp<volField> using inverse distance weighting
161  // returning pointField
162  template<class Type>
164  (
166  ) const;
167 
168  //- Interpolate volField using inverse distance weighting
169  // returning pointField with the same patchField types. Assigns
170  // to any fixedValue boundary conditions to make them consistent
171  // with internal field
172  template<class Type>
174  (
176  const wordList& patchFieldTypes
177  ) const;
178 
179  //- Interpolate tmp<volField> using inverse distance weighting
180  // returning pointField with the same patchField types. Assigns
181  // to any fixedValue boundary conditions to make them consistent
182  // with internal field
183  template<class Type>
185  (
187  const wordList& patchFieldTypes
188  ) const;
189 
190 
191  // Low level
192 
193  //- Interpolate internal field from volField to pointField
194  // using inverse distance weighting
195  template<class Type>
197  (
200  ) const;
201 
202  //- Interpolate boundary field without applying constraints/boundary
203  // conditions
204  template<class Type>
206  (
209  ) const;
210 
211  //- Interpolate boundary with constraints/boundary conditions
212  template<class Type>
214  (
217  const bool overrideFixedValue
218  ) const;
219 
220  //- Interpolate from volField to pointField
221  // using inverse distance weighting
222  template<class Type>
223  void interpolate
224  (
227  ) const;
228 
229  //- Interpolate volField using inverse distance weighting
230  // returning pointField with name. Optionally caches
231  template<class Type>
233  (
235  const word& name,
236  const bool cache
237  ) const;
238 
239 
240  // Interpolation for displacement (applies 2D correction)
241 
242  //- Interpolate from volField to pointField
243  // using inverse distance weighting
245  (
246  const volVectorField&,
248  ) const;
249 
250 };
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 } // End namespace Foam
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 #ifdef NoRepository
260  #include "volPointInterpolate.C"
261 #endif
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 #endif
266 
267 // ************************************************************************* //
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
const word & name() const
Return name.
Definition: IOobject.H:297
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
Interpolate boundary field without applying constraints/boundary.
Generic GeometricField class.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
ClassName("volPointInterpolation")
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:85
A class for handling words, derived from string.
Definition: word.H:59
void interpolateInternalField(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate internal field from volField to pointField.
Interpolate from cell centres to points (vertices) using inverse distance weighting.
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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
bool movePoints()
Correct weighting factors for moving mesh.
Namespace for OpenFOAM.