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-2024 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 "volFieldsFwd.H"
48 #include "pointFieldsFwd.H"
49 #include "scalarField.H"
50 #include "vectorField.H"
51 #include "Map.H"
52 #include "DynamicList.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 class face;
60 class polyMesh;
61 class pointMesh;
62 
63 /*---------------------------------------------------------------------------*\
64  Class pointMVCWeight Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 class pointMVCWeight
68 {
69 protected:
70 
71  // Protected data
72 
73  //- Cell index
74  const label cellIndex_;
75 
76  //- Weights applied to cell vertices
78 
79 
80  // Protected Member Functions
81 
82  //- Calculate weights from single face's vertices only
83  void calcWeights
84  (
85  const Map<label>& toLocal,
86  const face& f,
87  const DynamicList<point>& u,
88  const scalarField& dist,
90  ) const;
91 
92  //- Calculate weights from all cell's vertices
93  void calcWeights
94  (
95  const polyMesh& mesh,
96  const labelList& toGlobal,
97  const Map<label>& toLocal,
98  const vector& position,
99  const vectorField& uVec,
100  const scalarField& dist,
102  ) const;
103 
104 public:
105 
106  //- Type information
107  ClassName("pointMVCWeight");
108 
109  //- Tolerance used in calculating barycentric co-ordinates
110  // (applied to normalised values)
111  static scalar tol;
112 
113 
114  // Constructors
115 
116  //- Construct from components
118  (
119  const polyMesh& mesh,
120  const vector& position,
121  const label celli,
122  const label facei = -1
123  );
124 
125 
126  // Member Functions
127 
128  //- Cell index
129  inline label cell() const
130  {
131  return cellIndex_;
132  }
133 
134  //- Interpolation weights (in order of cellPoints)
135  inline const scalarField& weights() const
136  {
137  return weights_;
138  }
139 
140  //- Interpolate field
141  template<class Type>
142  inline Type interpolate
143  (
144  const PointField<Type>& psip
145  ) const;
146 };
147 
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 } // End namespace Foam
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 #include "pointMVCWeightI.H"
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 #endif
160 
161 // ************************************************************************* //
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
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)