volPointInterpolation.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  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::volPointInterpolation
29 
30 Description
31  Interpolate from cell centres to points (vertices) using inverse distance
32  weighting
33 
34 SourceFiles
35  volPointInterpolation.C
36  volPointInterpolate.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef volPointInterpolation_H
41 #define volPointInterpolation_H
42 
43 #include "MeshObject.H"
44 #include "scalarList.H"
45 #include "volFields.H"
46 #include "pointFields.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 // Forward Declarations
54 class fvMesh;
55 class pointMesh;
56 
57 /*---------------------------------------------------------------------------*\
58  Class volPointInterpolation Declaration
59 \*---------------------------------------------------------------------------*/
60 
62 :
63  public MeshObject<fvMesh, UpdateableMeshObject, volPointInterpolation>
64 {
65  // Private Typedefs
66 
67  typedef MeshObject
68  <
69  fvMesh,
73 
74 
75  // Private Data
76 
77  //- Interpolation scheme weighting factor array.
78  scalarListList pointWeights_;
79 
80 
81  // Boundary handling
82 
83  //- Boundary addressing
84  autoPtr<primitivePatch> boundaryPtr_;
85 
86  //- Per boundary face whether is on non-coupled, non-empty patch
87  bitSet boundaryIsPatchFace_;
88 
89  //- Per mesh(!) point whether is on non-coupled, non-empty patch (on
90  // any processor)
91  bitSet isPatchPoint_;
92 
93  //- Per boundary point the weights per pointFaces.
94  scalarListList boundaryPointWeights_;
95 
96  //- Whether mesh has any non-processor coupled patches
97  bool hasSeparated_;
98 
99  //- Optional scaling for coupled point patch fields using the
100  // swapAddSeparated to add additional contributions (e.g. cyclicAMI).
101  // (note: one-to-one matching (processor) uses a different mechanism)
102  tmp<scalarField> normalisationPtr_;
103 
104 
105  // Private Member Functions
106 
107  //- Constructor helper: check if any non-processor coupled patches
108  static bool hasSeparated(const pointMesh& pMesh);
109 
110  //- Construct addressing over all boundary faces
111  void calcBoundaryAddressing();
112 
113  //- Make weights for internal and coupled-only boundarypoints
114  void makeInternalWeights(scalarField& sumWeights);
115 
116  //- Make weights for points on uncoupled patches
117  void makeBoundaryWeights(scalarField& sumWeights);
118 
119  //- Construct all point weighting factors
120  void makeWeights();
121 
122  //- Helper: interpolate oneField onto boundary points only. Used to
123  // correct normalisation
124  void interpolateOne
125  (
126  const tmp<scalarField>& tnormalisation,
127  pointScalarField& pf
128  ) const;
129 
130  //- Helper: push master point data to collocated points
131  template<class Type>
132  void pushUntransformedData(List<Type>&) const;
133 
134  //- Get boundary field in same order as boundary faces. Field is
135  // zero on all coupled and empty patches
136  template<class Type>
137  tmp<Field<Type>> flatBoundaryField
138  (
140  ) const;
141 
142  //- Add separated contributions
143  template<class Type>
144  void addSeparated
145  (
147  ) const;
148 
149  //- No copy construct
151 
152  //- No copy assignment
153  void operator=(const volPointInterpolation&) = delete;
154 
155 
156 public:
157 
158  // Declare name of the class and its debug switch
159  ClassName("volPointInterpolation");
160 
161 
162  // Constructors
163 
164  //- Construct given fvMesh
165  explicit volPointInterpolation(const fvMesh&);
166 
167 
168  //- Destructor
170 
171 
172  // Member functions
173 
174  // Edit
175 
176  //- Update mesh topology using the morph engine
177  void updateMesh(const mapPolyMesh&);
178 
179  //- Correct weighting factors for moving mesh.
180  bool movePoints();
181 
182 
183  // Interpolation Functions
184 
185  //- Interpolate volField using inverse distance weighting
186  // \return pointField
187  template<class Type>
189  (
191  ) const;
192 
193  //- Interpolate tmp<volField> using inverse distance weighting
194  // \return pointField
195  template<class Type>
197  (
199  ) const;
200 
201  //- Interpolate volField using inverse distance weighting
202  // returning pointField with the same patchField types. Assigns
203  // to any fixedValue boundary conditions to make them consistent
204  // with internal field
205  template<class Type>
207  (
209  const wordList& patchFieldTypes
210  ) const;
211 
212  //- Interpolate tmp<volField> using inverse distance weighting
213  // returning pointField with the same patchField types. Assigns
214  // to any fixedValue boundary conditions to make them consistent
215  // with internal field
216  template<class Type>
218  (
220  const wordList& patchFieldTypes
221  ) const;
222 
223  //- Interpolate dimensionedField using inverse distance weighting
224  // \return pointField
225  template<class Type>
227  (
229  ) const;
230 
231  //- Interpolate tmp<dimensionedField> using inverse distance
232  // weighting returning pointField
233  template<class Type>
235  (
237  ) const;
238 
239 
240  // Low level
241 
242  //- Interpolate internal field from volField to pointField
243  // using inverse distance weighting
244  template<class Type>
246  (
249  ) const;
250 
251  //- Interpolate boundary field without applying constraints/boundary
252  // conditions
253  template<class Type>
255  (
257 
259  ) const;
260 
261  //- Interpolate boundary with constraints/boundary conditions
262  template<class Type>
264  (
267  const bool overrideFixedValue
268  ) const;
269 
270  //- Interpolate from volField to pointField
271  // using inverse distance weighting
272  template<class Type>
273  void interpolate
274  (
277  ) const;
278 
279  //- Interpolate volField using inverse distance weighting
280  // returning pointField with name. Optionally caches
281  template<class Type>
283  (
285  const word& name,
286  const bool cache
287  ) const;
288 
289  //- Interpolate dimensioned internal field from cells to points
290  // using inverse distance weighting
291  template<class Type>
293  (
296  ) const;
297 
298  //- Interpolate dimensionedField using inverse distance weighting
299  // returning pointField with name. Optionally caches
300  template<class Type>
302  (
304  const word& name,
305  const bool cache
306  ) const;
307 
308  // Interpolation for displacement (applies 2D correction)
309 
310  //- Interpolate from volField to pointField
311  // using inverse distance weighting
313  (
314  const volVectorField&,
316  ) const;
317 };
318 
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321 
322 } // End namespace Foam
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
326 #ifdef NoRepository
327  #include "volPointInterpolate.C"
328 #endif
329 
330 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331 
332 #endif
333 
334 // ************************************************************************* //
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.
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
Interpolate boundary field without applying constraints/boundary.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:205
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
ClassName("volPointInterpolation")
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:152
A class for handling words, derived from Foam::string.
Definition: word.H:63
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 interpolateDimensionedInternalField(const DimensionedField< Type, volMesh > &vf, DimensionedField< Type, pointMesh > &pf) const
Interpolate dimensioned internal field from cells to points.
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:58
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool movePoints()
Correct weighting factors for moving mesh.
Namespace for OpenFOAM.
UpdateableMeshObject(const word &objName, const objectRegistry &obr)
Construct from name and instance on registry.
Definition: MeshObject.H:354