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 | 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 \*---------------------------------------------------------------------------*/
27 
28 #include "pointMVCWeight.H"
29 #include "polyMesh.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineDebugSwitchWithName(pointMVCWeight, "pointMVCWeight", 0);
36 }
37 
38 Foam::scalar Foam::pointMVCWeight::tol(SMALL);
39 
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 (
45  const Map<label>& toLocal,
46  const face& f,
47  const DynamicList<point>& u,
48  const scalarField& dist,
49  scalarField& weights
50 ) const
51 {
52  weights.setSize(toLocal.size());
53  weights = 0.0;
54 
55  scalarField theta(f.size());
56 
57  // recompute theta, the theta computed previously are not robust
58  forAll(f, j)
59  {
60  label jPlus1 = f.fcIndex(j);
61  scalar l = mag(u[j] - u[jPlus1]);
62  theta[j] = 2.0*Foam::asin(l/2.0);
63  }
64 
65  scalar sumWeight = 0;
66  forAll(f, j)
67  {
68  label pid = toLocal[f[j]];
69  label jMin1 = f.rcIndex(j);
70  weights[pid] =
71  1.0
72  / dist[pid]
73  * (Foam::tan(theta[jMin1]/2.0) + Foam::tan(theta[j]/2.0));
74  sumWeight += weights[pid];
75  }
76 
77  if (sumWeight >= tol)
78  {
79  weights /= sumWeight;
80  }
81 }
82 
83 
85 (
86  const polyMesh& mesh,
87  const labelList& toGlobal,
88  const Map<label>& toLocal,
89  const vector& position,
90  const vectorField& uVec,
91  const scalarField& dist,
92  scalarField& weights
93 ) const
94 {
95  // Loop over all triangles of all polygons of cell to compute weights
97  DynamicList<scalar> theta(100);
98  DynamicList<point> u(100);
99 
100  const Foam::cell& cFaces = mesh.cells()[cellIndex_];
101 
102  forAll(cFaces, iter)
103  {
104  label facei = cFaces[iter];
105  const face& f = mesh.faces()[facei];
106 
107  //Pout<< "face:" << facei << " at:"
108  // << pointField(mesh.points(), f)
109  // << endl;
110 
111  // Collect the uVec for the face
112  forAll(f, j)
113  {
114  u(j) = uVec[toLocal[f[j]]];
115  }
116 
117  vector v(Zero);
118  forAll(f, j)
119  {
120  label jPlus1 = f.fcIndex(j);
121  //Pout<< " uj:" << u[j] << " ujPlus1:" << u[jPlus1] << endl;
122 
123  vector temp = u[j] ^ u[jPlus1];
124 
125  scalar magTemp = mag(temp);
126 
127  if (magTemp < VSMALL)
128  {
129  continue;
130  }
131 
132  temp /= magTemp;
133 
134  //Pout<< " uj:" << u[j] << " ujPlus1:" << u[jPlus1]
135  // << " temp:" << temp << endl;
136 
137  scalar l = min(mag(u[j] - u[jPlus1]), 2.0);
138  scalar angle = 2.0*Foam::asin(l/2.0);
139 
140  //Pout<< " j:" << j << " l:" << l << " angle:" << angle << endl;
141 
142  v += 0.5*angle*temp;
143  }
144 
145  scalar vNorm = mag(v);
146  v /= vNorm;
147 
148  if ((v & u[0]) < 0)
149  {
150  v = -v;
151  }
152 
153  //Pout<< " v:" << v << endl;
154 
155  // angles between edges
156  forAll(f, j)
157  {
158  label jPlus1 = f.fcIndex(j);
159  //Pout<< " uj:" << u[j] << " ujPlus1:" << u[jPlus1] << endl;
160 
161  const vector n0 = normalised(u[j] ^ v);
162  const vector n1 = normalised(u[jPlus1] ^ v);
163 
164  scalar l = min(mag(n0 - n1), 2.0);
165  //Pout<< " l:" << l << endl;
166  alpha(j) = 2.0*Foam::asin(l/2.0);
167 
168  vector temp = n0^n1;
169  if ((temp&v) < 0.0)
170  {
171  alpha[j] = -alpha[j];
172  }
173 
174  l = min(mag(u[j] - v), 2.0);
175  //Pout<< " l:" << l << endl;
176  theta(j) = 2.0*Foam::asin(l/2.0);
177  }
178 
179 
180  bool outlierFlag = false;
181  forAll(f, j)
182  {
183  if (mag(theta[j]) < tol)
184  {
185  outlierFlag = true;
186 
187  label pid = toLocal[f[j]];
188  weights[pid] += vNorm / dist[pid];
189  break;
190  }
191  }
192 
193  if (outlierFlag)
194  {
195  continue;
196  }
197 
198  scalar sum = 0.0;
199  forAll(f, j)
200  {
201  label jMin1 = f.rcIndex(j);
202  sum +=
203  1.0
204  / Foam::tan(theta[j])
205  * (Foam::tan(alpha[j]/2.0) + Foam::tan(alpha[jMin1]/2.0));
206  }
207 
208  // The special case when x lies on the polygon, handle it using 2D mvc.
209  // In the 2D case, alpha = theta
210  if (mag(sum) < tol)
211  {
212  // Calculate weights using face vertices only
213  calcWeights(toLocal, f, u, dist, weights);
214  return;
215  }
216 
217 
218  // Normal 3D case
219  forAll(f, j)
220  {
221  label pid = toLocal[f[j]];
222  label jMin1 = f.rcIndex(j);
223  weights[pid] +=
224  vNorm
225  / sum
226  / dist[pid]
227  / Foam::sin(theta[j])
228  * (Foam::tan(alpha[j]/2.0) + Foam::tan(alpha[jMin1]/2.0));
229  }
230  }
231 
232  // normalise weights
233  scalar sumWeight = sum(weights);
234 
235  if (mag(sumWeight) < tol)
236  {
237  return;
238  }
239  weights /= sumWeight;
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
244 
246 (
247  const polyMesh& mesh,
248  const vector& position,
249  const label cellIndex,
250  const label faceIndex
251 )
252 :
253  cellIndex_((cellIndex != -1) ? cellIndex : mesh.faceOwner()[faceIndex])
254 {
255  // Addressing - face vertices to local points and vice versa
256  const labelList& toGlobal = mesh.cellPoints()[cellIndex_];
257  Map<label> toLocal(2*toGlobal.size());
258  forAll(toGlobal, i)
259  {
260  toLocal.insert(toGlobal[i], i);
261  }
262 
263 
264  // Initialise weights
265  weights_.setSize(toGlobal.size());
266  weights_ = 0.0;
267 
268 
269  // Point-to-vertex vectors and distances
270  vectorField uVec(toGlobal.size());
271  scalarField dist(toGlobal.size());
272 
273  forAll(toGlobal, pid)
274  {
275  const point& pt = mesh.points()[toGlobal[pid]];
276 
277  uVec[pid] = pt-position;
278  dist[pid] = mag(uVec[pid]);
279 
280  // Special case: point is close to vertex
281  if (dist[pid] < tol)
282  {
283  weights_[pid] = 1.0;
284  return;
285  }
286  }
287 
288  // Project onto unit sphere
289  uVec /= dist;
290 
291 
292  if (faceIndex < 0)
293  {
294  // Face data not supplied
295  calcWeights
296  (
297  mesh,
298  toGlobal,
299  toLocal,
300  position,
301  uVec,
302  dist,
303 
304  weights_
305  );
306  }
307  else
308  {
309  DynamicList<point> u(100);
310  const face& f = mesh.faces()[faceIndex];
311  // Collect the uVec for the face
312  forAll(f, j)
313  {
314  u(j) = uVec[toLocal[f[j]]];
315  }
316 
317  // Calculate weights for face only
318  calcWeights(toLocal, f, u, dist, weights_);
319  }
320 }
321 
322 
323 // ************************************************************************* //
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::defineDebugSwitchWithName
defineDebugSwitchWithName(pointMVCWeight, "pointMVCWeight", 0)
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< point >
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::Map< label >
polyMesh.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
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
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Field< scalar >
Foam::pid
pid_t pid()
Return the PID of this process.
Definition: MSwindows.C:330
Foam::pointMVCWeight::pointMVCWeight
pointMVCWeight(const polyMesh &mesh, const vector &position, const label celli, const label facei=-1)
Construct from components.
Definition: pointMVCWeight.C:246
pointMVCWeight.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267