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 -------------------------------------------------------------------------------
10 License
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 
26 Class
27  Foam::pointMVCWeight
28 
29 Description
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 
38 
39 
40 SourceFiles
41  pointMVCWeight.C
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #ifndef pointMVCWeight_H
46 #define pointMVCWeight_H
47 
48 #include "scalarField.H"
49 #include "vectorField.H"
50 #include "Map.H"
51 #include "DynamicList.H"
52 #include "point.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 class polyMesh;
60 class pointMesh;
61 template<class T> class pointPatchField;
62 template<class Type, template<class> class PatchField, class GeoMesh>
63 class GeometricField;
64 class face;
65 
66 /*---------------------------------------------------------------------------*\
67  Class pointMVCWeight Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class pointMVCWeight
71 {
72 protected:
73 
74  // Protected data
75 
76  //- Cell index
77  const label cellIndex_;
78 
79  //- Weights applied to cell vertices
81 
82 
83  // Protected Member Functions
84 
85  //- Calculate weights from single face's vertices only
86  void calcWeights
87  (
88  const Map<label>& toLocal,
89  const face& f,
90  const DynamicList<point>& u,
91  const scalarField& dist,
93  ) const;
94 
95  //- Calculate weights from all cell's vertices
96  void calcWeights
97  (
98  const polyMesh& mesh,
99  const labelList& toGlobal,
100  const Map<label>& toLocal,
101  const vector& position,
102  const vectorField& uVec,
103  const scalarField& dist,
105  ) const;
106 
107 public:
108 
109  //- Debug switch
110  static int debug;
111 
112  //- Tolerance used in calculating barycentric coordinates
113  // (applied to normalised values)
114  static scalar tol;
115 
116 
117  // Constructors
118 
119  //- Construct from components
121  (
122  const polyMesh& mesh,
123  const vector& position,
124  const label celli,
125  const label facei = -1
126  );
127 
128 
129  // Member Functions
130 
131  //- Cell index
132  inline label cell() const
133  {
134  return cellIndex_;
135  }
136 
137  //- Interpolation weights (in order of cellPoints)
138  inline const scalarField& weights() const
139  {
140  return weights_;
141  }
142 
143  //- Interpolate field
144  template<class Type>
145  inline Type interpolate
146  (
148  ) const;
149 };
150 
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 } // End namespace Foam
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 #include "pointMVCWeightI.H"
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 #endif
163 
164 // ************************************************************************* //
Foam::pointMVCWeight::debug
static int debug
Debug switch.
Definition: pointMVCWeight.H:109
Foam::pointMVCWeight::weights_
scalarField weights_
Weights applied to cell vertices.
Definition: pointMVCWeight.H:79
scalarField.H
pointMVCWeightI.H
Foam::DynamicList< point >
point.H
Foam::Map< label >
Foam::pointMVCWeight::interpolate
Type interpolate(const GeometricField< Type, pointPatchField, pointMesh > &psip) const
Interpolate field.
Definition: pointMVCWeightI.H:34
Foam::pointPatchField
Abstract base class for point-mesh patch fields.
Definition: pointMVCWeight.H:60
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::pointMVCWeight::tol
static scalar tol
Tolerance used in calculating barycentric coordinates.
Definition: pointMVCWeight.H:113
Map.H
Foam::Field< scalar >
Foam::pointMVCWeight::pointMVCWeight
pointMVCWeight(const polyMesh &mesh, const vector &position, const label celli, const label facei=-1)
Construct from components.
Definition: pointMVCWeight.C:246
Foam::pointMVCWeight
Container to calculate weights for interpolating directly from vertices of cell using Mean Value Coor...
Definition: pointMVCWeight.H:69
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:48
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:51
Foam::pointMVCWeight::cellIndex_
const label cellIndex_
Cell index.
Definition: pointMVCWeight.H:76
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
vectorField.H
Foam::pointMVCWeight::calcWeights
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.
Definition: pointMVCWeight.C:44
Foam::pointMVCWeight::cell
label cell() const
Cell index.
Definition: pointMVCWeight.H:131
Foam::pointMVCWeight::weights
const scalarField & weights() const
Interpolation weights (in order of cellPoints)
Definition: pointMVCWeight.H:137
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
DynamicList.H
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53