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-2022 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 "point.H"
47 #include "pointFieldsFwd.H"
48 #include "scalarField.H"
49 #include "vectorField.H"
50 #include "Map.H"
51 #include "DynamicList.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 class face;
59 class polyMesh;
60 
61 /*---------------------------------------------------------------------------*\
62  Class pointMVCWeight Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 class pointMVCWeight
66 {
67 protected:
68 
69  // Protected data
70 
71  //- Cell index
72  const label cellIndex_;
73 
74  //- Weights applied to cell vertices
76 
77 
78  // Protected Member Functions
79 
80  //- Calculate weights from single face's vertices only
81  void calcWeights
82  (
83  const Map<label>& toLocal,
84  const face& f,
85  const DynamicList<point>& u,
86  const scalarField& dist,
88  ) const;
89 
90  //- Calculate weights from all cell's vertices
91  void calcWeights
92  (
93  const polyMesh& mesh,
94  const labelList& toGlobal,
95  const Map<label>& toLocal,
96  const vector& position,
97  const vectorField& uVec,
98  const scalarField& dist,
100  ) const;
101 
102 public:
103 
104  //- Type information
105  ClassName("pointMVCWeight");
106 
107  //- Tolerance used in calculating barycentric co-ordinates
108  // (applied to normalised values)
109  static scalar tol;
110 
111 
112  // Constructors
113 
114  //- Construct from components
116  (
117  const polyMesh& mesh,
118  const vector& position,
119  const label celli,
120  const label facei = -1
121  );
122 
123 
124  // Member Functions
125 
126  //- Cell index
127  inline label cell() const
128  {
129  return cellIndex_;
130  }
131 
132  //- Interpolation weights (in order of cellPoints)
133  inline const scalarField& weights() const
134  {
135  return weights_;
136  }
137 
138  //- Interpolate field
139  template<class Type>
140  inline Type interpolate
141  (
142  const PointField<Type>& psip
143  ) const;
144 };
145 
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 } // End namespace Foam
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 #include "pointMVCWeightI.H"
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 #endif
158 
159 // ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coor...
const label cellIndex_
Cell index.
static scalar tol
Tolerance used in calculating barycentric co-ordinates.
label cell() const
Cell index.
scalarField weights_
Weights applied to cell vertices.
const scalarField & weights() const
Interpolation weights (in order of cellPoints)
Type interpolate(const PointField< Type > &psip) const
Interpolate field.
pointMVCWeight(const polyMesh &mesh, const vector &position, const label celli, const label facei=-1)
Construct from components.
void calcWeights(const Map< label > &toLocal, const face &f, const DynamicList< point > &u, const scalarField &dist, scalarField &weights) const
Calculate weights from single face's vertices only.
ClassName("pointMVCWeight")
Type information.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Namespace for OpenFOAM.
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
labelList f(nPoints)