pointMVCWeight.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-2020 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::pointMVCWeight
26 
27 Description
28  Container to calculate weights for interpolating directly from vertices
29  of cell using Mean Value Coordinates.
30 
31  Based on (VTK's vtkMeanValueCoordinatesInterpolator's) implementation
32  of "Spherical Barycentric Coordinates"
33  2006 paper Eurographics Symposium on Geometry Processing
34  by Torsten Langer, Alexander Belyaev and Hans-Peter Seide
35 
36 
37 
38 SourceFiles
39  pointMVCWeight.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef pointMVCWeight_H
44 #define pointMVCWeight_H
45 
46 #include "scalarField.H"
47 #include "vectorField.H"
48 #include "Map.H"
49 #include "DynamicList.H"
50 #include "point.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 class polyMesh;
58 class pointMesh;
59 template<class T> class pointPatchField;
60 template<class Type, template<class> class PatchField, class GeoMesh>
61 class GeometricField;
62 class face;
63 
64 /*---------------------------------------------------------------------------*\
65  Class pointMVCWeight Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class pointMVCWeight
69 {
70 protected:
71 
72  // Protected data
73 
74  //- Cell index
75  const label cellIndex_;
76 
77  //- Weights applied to cell vertices
79 
80 
81  // Protected Member Functions
82 
83  //- Calculate weights from single face's vertices only
84  void calcWeights
85  (
86  const Map<label>& toLocal,
87  const face& f,
88  const DynamicList<point>& u,
89  const scalarField& dist,
91  ) const;
92 
93  //- Calculate weights from all cell's vertices
94  void calcWeights
95  (
96  const polyMesh& mesh,
97  const labelList& toGlobal,
98  const Map<label>& toLocal,
99  const vector& position,
100  const vectorField& uVec,
101  const scalarField& dist,
103  ) const;
104 
105 public:
106 
107  //- Type information
108  ClassName("pointMVCWeight");
109 
110  //- Tolerance used in calculating barycentric co-ordinates
111  // (applied to normalised values)
112  static scalar tol;
113 
114 
115  // Constructors
116 
117  //- Construct from components
119  (
120  const polyMesh& mesh,
121  const vector& position,
122  const label celli,
123  const label facei = -1
124  );
125 
126 
127  // Member Functions
128 
129  //- Cell index
130  inline label cell() const
131  {
132  return cellIndex_;
133  }
134 
135  //- Interpolation weights (in order of cellPoints)
136  inline const scalarField& weights() const
137  {
138  return weights_;
139  }
140 
141  //- Interpolate field
142  template<class Type>
143  inline Type interpolate
144  (
146  ) const;
147 };
148 
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 } // End namespace Foam
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 #include "pointMVCWeightI.H"
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 #endif
161 
162 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
ClassName("pointMVCWeight")
Type information.
Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coor...
Generic GeometricField class.
static scalar tol
Tolerance used in calculating barycentric co-ordinates.
scalarField weights_
Weights applied to cell vertices.
Type interpolate(const GeometricField< Type, pointPatchField, pointMesh > &psip) const
Interpolate field.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Abstract base class for point-mesh patch fields.
dynamicFvMesh & mesh
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
pointMVCWeight(const polyMesh &mesh, const vector &position, const label celli, const label facei=-1)
Construct from components.
const label cellIndex_
Cell index.
void calcWeights(const Map< label > &toLocal, const face &f, const DynamicList< point > &u, const scalarField &dist, scalarField &weights) const
Calculate weights from single face&#39;s vertices only.
labelList f(nPoints)
label cell() const
Cell index.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
const scalarField & weights() const
Interpolation weights (in order of cellPoints)
Namespace for OpenFOAM.