pointMVCWeight.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2012 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 \*---------------------------------------------------------------------------*/
25 
26 #include "pointMVCWeight.H"
27 #include "polyMesh.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 defineDebugSwitchWithName(pointMVCWeight, "pointMVCWeight", 0);
34 }
35 
36 Foam::scalar Foam::pointMVCWeight::tol(SMALL);
37 
38 
39 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
40 
42 (
43  const Map<label>& toLocal,
44  const face& f,
45  const DynamicList<point>& u,
46  const scalarField& dist,
47  scalarField& weights
48 ) const
49 {
50  weights.setSize(toLocal.size());
51  weights = 0.0;
52 
53  scalarField theta(f.size());
54 
55  // recompute theta, the theta computed previously are not robust
56  forAll(f, j)
57  {
58  label jPlus1 = f.fcIndex(j);
59  scalar l = mag(u[j] - u[jPlus1]);
60  theta[j] = 2.0*Foam::asin(l/2.0);
61  }
62 
63  scalar sumWeight = 0;
64  forAll(f, j)
65  {
66  label pid = toLocal[f[j]];
67  label jMin1 = f.rcIndex(j);
68  weights[pid] =
69  1.0
70  / dist[pid]
71  * (Foam::tan(theta[jMin1]/2.0) + Foam::tan(theta[j]/2.0));
72  sumWeight += weights[pid];
73  }
74 
75  if (sumWeight >= tol)
76  {
77  weights /= sumWeight;
78  }
79 }
80 
81 
83 (
84  const polyMesh& mesh,
85  const labelList& toGlobal,
86  const Map<label>& toLocal,
87  const vector& position,
88  const vectorField& uVec,
89  const scalarField& dist,
90  scalarField& weights
91 ) const
92 {
93  // Loop over all triangles of all polygons of cell to compute weights
95  DynamicList<scalar> theta(100);
96  DynamicList<point> u(100);
97 
98  const Foam::cell& cFaces = mesh.cells()[cellIndex_];
99 
100  forAll(cFaces, iter)
101  {
102  label faceI = cFaces[iter];
103  const face& f = mesh.faces()[faceI];
104 
105  //Pout<< "face:" << faceI << " at:"
106  // << pointField(mesh.points(), f)
107  // << endl;
108 
109  // Collect the uVec for the face
110  forAll(f, j)
111  {
112  u(j) = uVec[toLocal[f[j]]];
113  }
114 
115  vector v(point::zero);
116  forAll(f, j)
117  {
118  label jPlus1 = f.fcIndex(j);
119  //Pout<< " uj:" << u[j] << " ujPlus1:" << u[jPlus1] << endl;
120 
121  vector temp = u[j] ^ u[jPlus1];
122 
123  scalar magTemp = mag(temp);
124 
125  if (magTemp < VSMALL)
126  {
127  continue;
128  }
129 
130  temp /= magTemp;
131 
132  //Pout<< " uj:" << u[j] << " ujPlus1:" << u[jPlus1]
133  // << " temp:" << temp << endl;
134 
135  scalar l = min(mag(u[j] - u[jPlus1]), 2.0);
136  scalar angle = 2.0*Foam::asin(l/2.0);
137 
138  //Pout<< " j:" << j << " l:" << l << " angle:" << angle << endl;
139 
140  v += 0.5*angle*temp;
141  }
142 
143  scalar vNorm = mag(v);
144  v /= vNorm;
145 
146  // Make sure v points towards the polygon
147  //if (((v&u[0]) < 0) != (mesh.faceOwner()[faceI] != cellIndex_))
148  //{
149  // FatalErrorIn("pointMVCWeight::calcWeights(..)")
150  // << "v:" << v << " u[0]:" << u[0]
151  // << exit(FatalError);
152  //}
153 
154  if ((v & u[0]) < 0)
155  {
156  v = -v;
157  }
158 
159  //Pout<< " v:" << v << endl;
160 
161  // angles between edges
162  forAll(f, j)
163  {
164  label jPlus1 = f.fcIndex(j);
165  //Pout<< " uj:" << u[j] << " ujPlus1:" << u[jPlus1] << endl;
166 
167  vector n0 = u[j]^v;
168  n0 /= mag(n0);
169  vector n1 = u[jPlus1]^v;
170  n1 /= mag(n1);
171 
172  scalar l = min(mag(n0 - n1), 2.0);
173  //Pout<< " l:" << l << endl;
174  alpha(j) = 2.0*Foam::asin(l/2.0);
175 
176  vector temp = n0^n1;
177  if ((temp&v) < 0.0)
178  {
179  alpha[j] = -alpha[j];
180  }
181 
182  l = min(mag(u[j] - v), 2.0);
183  //Pout<< " l:" << l << endl;
184  theta(j) = 2.0*Foam::asin(l/2.0);
185  }
186 
187 
188  bool outlierFlag = false;
189  forAll(f, j)
190  {
191  if (mag(theta[j]) < tol)
192  {
193  outlierFlag = true;
194 
195  label pid = toLocal[f[j]];
196  weights[pid] += vNorm / dist[pid];
197  break;
198  }
199  }
200 
201  if (outlierFlag)
202  {
203  continue;
204  }
205 
206  scalar sum = 0.0;
207  forAll(f, j)
208  {
209  label jMin1 = f.rcIndex(j);
210  sum +=
211  1.0
212  / Foam::tan(theta[j])
213  * (Foam::tan(alpha[j]/2.0) + Foam::tan(alpha[jMin1]/2.0));
214  }
215 
216  // The special case when x lies on the polygon, handle it using 2D mvc.
217  // In the 2D case, alpha = theta
218  if (mag(sum) < tol)
219  {
220  // Calculate weights using face vertices only
221  calcWeights(toLocal, f, u, dist, weights);
222  return;
223  }
224 
225 
226  // Normal 3D case
227  forAll(f, j)
228  {
229  label pid = toLocal[f[j]];
230  label jMin1 = f.rcIndex(j);
231  weights[pid] +=
232  vNorm
233  / sum
234  / dist[pid]
235  / Foam::sin(theta[j])
236  * (Foam::tan(alpha[j]/2.0) + Foam::tan(alpha[jMin1]/2.0));
237  }
238  }
239 
240  // normalise weights
241  scalar sumWeight = sum(weights);
242 
243  if (mag(sumWeight) < tol)
244  {
245  return;
246  }
247  weights /= sumWeight;
248 }
249 
250 
251 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
252 
254 (
255  const polyMesh& mesh,
256  const vector& position,
257  const label cellIndex,
258  const label faceIndex
259 )
260 :
261  cellIndex_((cellIndex != -1) ? cellIndex : mesh.faceOwner()[faceIndex])
262 {
263  // Addressing - face vertices to local points and vice versa
264  const labelList& toGlobal = mesh.cellPoints()[cellIndex_];
265  Map<label> toLocal(2*toGlobal.size());
266  forAll(toGlobal, i)
267  {
268  toLocal.insert(toGlobal[i], i);
269  }
270 
271 
272  // Initialise weights
273  weights_.setSize(toGlobal.size());
274  weights_ = 0.0;
275 
276 
277  // Point-to-vertex vectors and distances
278  vectorField uVec(toGlobal.size());
279  scalarField dist(toGlobal.size());
280 
281  forAll(toGlobal, pid)
282  {
283  const point& pt = mesh.points()[toGlobal[pid]];
284 
285  uVec[pid] = pt-position;
286  dist[pid] = mag(uVec[pid]);
287 
288  // Special case: point is close to vertex
289  if (dist[pid] < tol)
290  {
291  weights_[pid] = 1.0;
292  return;
293  }
294  }
295 
296  // Project onto unit sphere
297  uVec /= dist;
298 
299 
300  if (faceIndex < 0)
301  {
302  // Face data not supplied
303  calcWeights
304  (
305  mesh,
306  toGlobal,
307  toLocal,
308  position,
309  uVec,
310  dist,
311 
312  weights_
313  );
314  }
315  else
316  {
317  DynamicList<point> u(100);
318  const face& f = mesh.faces()[faceIndex];
319  // Collect the uVec for the face
320  forAll(f, j)
321  {
322  u(j) = uVec[toLocal[f[j]]];
323  }
324 
325  // Calculate weights for face only
326  calcWeights(toLocal, f, u, dist, weights_);
327  }
328 }
329 
330 
331 // ************************************************************************* //
const labelListList & cellPoints() const
dimensioned< scalar > mag(const dimensioned< Type > &)
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
defineDebugSwitchWithName(pointMVCWeight,"pointMVCWeight", 0)
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
const cellList & cells() const
Namespace for OpenFOAM.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
void setSize(const label)
Reset size of List.
Definition: List.C:318
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1035
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
#define forAll(list, i)
Definition: UList.H:421
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1073
dimensionedScalar tan(const dimensionedScalar &ds)
pid_t pid()
Return the PID of this process.
Definition: POSIX.C:78
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
static const Vector zero
Definition: Vector.H:80
static scalar tol
Tolerance used in calculating barycentric co-ordinates.
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1060
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
pointMVCWeight(const polyMesh &mesh, const vector &position, const label cellI, const label faceI=-1)
Construct from components.
dimensionedScalar asin(const dimensionedScalar &ds)
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.
dimensionedScalar sin(const dimensionedScalar &ds)