pointMVCWeight.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::pointMVCWeight
28
29Description
30 Container to calculate weights for interpolating directly from vertices
31 of cell using Mean Value Coordinates.
32
33 Based on (VTK's vtkMeanValueCoordinatesInterpolator's) implementation
34 of "Spherical Barycentric Coordinates"
35 2006 paper Eurographics Symposium on Geometry Processing
36 by Torsten Langer, Alexander Belyaev and Hans-Peter Seide
37
38SourceFiles
39 pointMVCWeight.C
40
41\*---------------------------------------------------------------------------*/
42
43#ifndef Foam_pointMVCWeight_H
44#define Foam_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
54namespace Foam
55{
56
57// Forward Declarations
58class face;
59class pointMesh;
60class polyMesh;
61template<class Type> class pointPatchField;
62template<class Type, template<class> class PatchField, class GeoMesh>
63class GeometricField;
64
65/*---------------------------------------------------------------------------*\
66 Class pointMVCWeight Declaration
67\*---------------------------------------------------------------------------*/
70{
71protected:
72
73 // Protected Data
74
75 //- Cell index
76 const label cellIndex_;
77
78 //- Weights applied to cell vertices
80
81
82 // Protected Member Functions
83
84 //- Calculate weights from single face's vertices only
85 void calcWeights
86 (
87 const Map<label>& toLocal,
88 const face& f,
89 const DynamicList<point>& u,
90 const scalarField& dist,
92 ) const;
93
94 //- Calculate weights from all cell's vertices
95 void calcWeights
96 (
97 const polyMesh& mesh,
98 const labelList& toGlobal,
99 const Map<label>& toLocal,
100 const vector& position,
101 const vectorField& uVec,
102 const scalarField& dist,
104 ) const;
105
106public:
107
108 //- Debug switch
109 static int debug;
110
111 //- Tolerance used in calculating barycentric coordinates
112 // (applied to normalised values)
113 static scalar tol;
114
115
116 // Constructors
117
118 //- Construct from components
120 (
121 const polyMesh& mesh,
122 const vector& position,
123 const label celli,
124 const label facei = -1
125 );
126
127
128 // Member Functions
129
130 //- Cell index
131 label cell() const noexcept
132 {
133 return cellIndex_;
134 }
135
136 //- Interpolation weights (in order of cellPoints)
137 const scalarField& weights() const noexcept
138 {
139 return weights_;
140 }
141
142 //- Interpolate field
143 template<class Type>
144 inline Type interpolate
145 (
147 ) const;
148};
149
150
151// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152
153} // End namespace Foam
154
155// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156
157#include "pointMVCWeightI.H"
158
159// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160
161#endif
162
163// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
Generic GeometricField class.
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coor...
Type interpolate(const GeometricField< Type, pointPatchField, pointMesh > &psip) const
Interpolate field.
const label cellIndex_
Cell index.
static scalar tol
Tolerance used in calculating barycentric coordinates.
scalarField weights_
Weights applied to cell vertices.
label cell() const noexcept
Cell index.
static int debug
Debug switch.
const scalarField & weights() const noexcept
Interpolation weights (in order of cellPoints)
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.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
dynamicFvMesh & mesh
Namespace for OpenFOAM.
const direction noexcept
Definition: Scalar.H:223
labelList f(nPoints)