All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2019 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 
117 public:
118 
119  // Declare name of the class and its debug switch
120  ClassName("volPointInterpolation");
121 
122 
123  // Constructors
124 
125  //- Constructor given fvMesh and pointMesh.
126  explicit volPointInterpolation(const fvMesh&);
127 
128  //- Disallow default bitwise copy construction
130 
131 
132  //- Destructor
134 
135 
136  // Member Functions
137 
138  // Edit
139 
140  //- Update mesh topology using the morph engine
141  void updateMesh(const mapPolyMesh&);
142 
143  //- Correct weighting factors for moving mesh.
144  bool movePoints();
145 
146 
147  // Interpolation functions
148 
149  //- Interpolate volField using inverse distance weighting
150  // returning pointField
151  template<class Type>
153  (
155  ) const;
156 
157  //- Interpolate tmp<volField> using inverse distance weighting
158  // returning pointField
159  template<class Type>
161  (
163  ) const;
164 
165  //- Interpolate volField using inverse distance weighting
166  // returning pointField with the same patchField types. Assigns
167  // to any fixedValue boundary conditions to make them consistent
168  // with internal field
169  template<class Type>
171  (
173  const wordList& patchFieldTypes
174  ) const;
175 
176  //- Interpolate tmp<volField> using inverse distance weighting
177  // returning pointField with the same patchField types. Assigns
178  // to any fixedValue boundary conditions to make them consistent
179  // with internal field
180  template<class Type>
182  (
184  const wordList& patchFieldTypes
185  ) const;
186 
187 
188  // Low level
189 
190  //- Interpolate internal field from volField to pointField
191  // using inverse distance weighting
192  template<class Type>
194  (
197  ) const;
198 
199  //- Interpolate boundary field without applying constraints/boundary
200  // conditions
201  template<class Type>
203  (
206  ) const;
207 
208  //- Interpolate boundary with constraints/boundary conditions
209  template<class Type>
211  (
214  const bool overrideFixedValue
215  ) const;
216 
217  //- Interpolate from volField to pointField
218  // using inverse distance weighting
219  template<class Type>
220  void interpolate
221  (
224  ) const;
225 
226  //- Interpolate volField using inverse distance weighting
227  // returning pointField with name. Optionally caches
228  template<class Type>
230  (
232  const word& name,
233  const bool cache
234  ) const;
235 
236 
237  // Interpolation for displacement (applies 2D correction)
238 
239  //- Interpolate from volField to pointField
240  // using inverse distance weighting
242  (
243  const volVectorField&,
245  ) const;
246 
247 
248  // Member Operators
249 
250  //- Disallow default bitwise assignment
251  void operator=(const volPointInterpolation&) = delete;
252 };
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 } // End namespace Foam
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 #ifdef NoRepository
262  #include "volPointInterpolate.C"
263 #endif
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 #endif
268 
269 // ************************************************************************* //
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:303
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:86
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.
volPointInterpolation(const fvMesh &)
Constructor given fvMesh and pointMesh.
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.
void operator=(const volPointInterpolation &)=delete
Disallow default bitwise assignment.