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 | Copyright (C) 2011-2014 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.
65 
66 
67  // Boundary handling
68 
69  //- Boundary addressing
71 
72  //- Per boundary face whether is on non-coupled patch
74 
75  //- Per mesh(!) point whether is on non-coupled patch (on any
76  // processor)
78 
79  //- Per boundary point the weights per pointFaces.
81 
82 
83  // Private Member Functions
84 
85  //- Construct addressing over all boundary faces
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>
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  //- Interpolate dimensionedField using inverse distance weighting
191  // returning pointField
192  template<class Type>
194  (
196  ) const;
197 
198  //- Interpolate tmp<dimensionedField> using inverse distance
199  // weighting returning pointField
200  template<class Type>
202  (
204  ) const;
205 
206 
207  // Low level
208 
209  //- Interpolate internal field from volField to pointField
210  // using inverse distance weighting
211  template<class Type>
213  (
216  ) const;
217 
218  //- Interpolate boundary field without applying constraints/boundary
219  // conditions
220  template<class Type>
222  (
225  ) const;
226 
227  //- Interpolate boundary with constraints/boundary conditions
228  template<class Type>
230  (
233  const bool overrideFixedValue
234  ) const;
235 
236  //- Interpolate from volField to pointField
237  // using inverse distance weighting
238  template<class Type>
239  void interpolate
240  (
243  ) const;
244 
245  //- Interpolate volField using inverse distance weighting
246  // returning pointField with name. Optionally caches
247  template<class Type>
249  (
251  const word& name,
252  const bool cache
253  ) const;
254 
255  //- Interpolate dimensioned internal field from cells to points
256  // using inverse distance weighting
257  template<class Type>
259  (
262  ) const;
263 
264  //- Interpolate dimensionedField using inverse distance weighting
265  // returning pointField with name. Optionally caches
266  template<class Type>
268  (
270  const word& name,
271  const bool cache
272  ) const;
273 
274  // Interpolation for displacement (applies 2D correction)
275 
276  //- Interpolate from volField to pointField
277  // using inverse distance weighting
279  (
280  const volVectorField&,
282  ) const;
283 
284 };
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace Foam
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 #ifdef NoRepository
294 # include "volPointInterpolate.C"
295 #endif
296 
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298 
299 #endif
300 
301 // ************************************************************************* //
volFields.H
Foam::volPointInterpolation::flatBoundaryField
tmp< Field< Type > > flatBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Get boundary field in same order as boundary faces. Field is.
Definition: volPointInterpolate.C:224
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::volPointInterpolation::ClassName
ClassName("volPointInterpolation")
Foam::volPointInterpolation::boundaryIsPatchFace_
boolList boundaryIsPatchFace_
Per boundary face whether is on non-coupled patch.
Definition: volPointInterpolation.H:72
Foam::volPointInterpolation::makeBoundaryWeights
void makeBoundaryWeights(scalarField &sumWeights)
Make weights for points on uncoupled patches.
Definition: volPointInterpolation.C:197
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::volPointInterpolation::volPointInterpolation
volPointInterpolation(const volPointInterpolation &)
Disallow default bitwise copy construct.
Foam::volPointInterpolation::boundaryPtr_
autoPtr< primitivePatch > boundaryPtr_
Boundary addressing.
Definition: volPointInterpolation.H:69
Foam::volPointInterpolation::updateMesh
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
Definition: volPointInterpolation.C:367
scalarList.H
Foam::volPointInterpolation::makeWeights
void makeWeights()
Construct all point weighting factors.
Definition: volPointInterpolation.C:245
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::volPointInterpolation::addSeparated
void addSeparated(GeometricField< Type, pointPatchField, pointMesh > &) const
Add separated contributions.
Definition: volPointInterpolate.C:85
Foam::volPointInterpolation::calcBoundaryAddressing
void calcBoundaryAddressing()
Construct addressing over all boundary faces.
Definition: volPointInterpolation.C:46
Foam::volPointInterpolation::operator=
void operator=(const volPointInterpolation &)
Disallow default bitwise assignment.
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::volPointInterpolation::~volPointInterpolation
~volPointInterpolation()
Destructor.
Definition: volPointInterpolation.C:361
Foam::volPointInterpolation::interpolate
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Definition: volPointInterpolate.C:502
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
volPointInterpolate.C
Foam::List< scalarList >
Foam::volPointInterpolation::movePoints
bool movePoints()
Correct weighting factors for moving mesh.
Definition: volPointInterpolation.C:373
Foam::volPointInterpolation::interpolateInternalField
void interpolateInternalField(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate internal field from volField to pointField.
Definition: volPointInterpolate.C:127
Foam::volPointInterpolation::interpolateDisplacement
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
Definition: volPointInterpolation.C:382
Foam::volPointInterpolation::makeInternalWeights
void makeInternalWeights(scalarField &sumWeights)
Make weights for internal and coupled-only boundarypoints.
Definition: volPointInterpolation.C:157
MeshObject.H
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::MeshObject
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:81
Foam::volPointInterpolation::boundaryPointWeights_
scalarListList boundaryPointWeights_
Per boundary point the weights per pointFaces.
Definition: volPointInterpolation.H:79
Foam::volPointInterpolation::isPatchPoint_
boolList isPatchPoint_
Per mesh(!) point whether is on non-coupled patch (on any.
Definition: volPointInterpolation.H:76
Foam::volPointInterpolation::pointWeights_
scalarListList pointWeights_
Interpolation scheme weighting factor array.
Definition: volPointInterpolation.H:63
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::volPointInterpolation::interpolateDimensionedInternalField
void interpolateDimensionedInternalField(const DimensionedField< Type, volMesh > &vf, DimensionedField< Type, pointMesh > &pf) const
Interpolate dimensioned internal field from cells to points.
Definition: volPointInterpolate.C:164
Foam::volPointInterpolation::interpolateBoundaryField
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
Interpolate boundary field without applying constraints/boundary.
Definition: volPointInterpolate.C:271
Foam::volPointInterpolation
Interpolate from cell centres to points (vertices) using inverse distance weighting.
Definition: volPointInterpolation.H:56
Foam::volPointInterpolation::pushUntransformedData
void pushUntransformedData(List< Type > &) const
Helper: push master point data to collocated points.
Definition: volPointInterpolate.C:42
pointFields.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51