primitiveMeshTools.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "primitiveMeshTools.H"
30 #include "syncTools.H"
31 #include "pyramidPointFaceRef.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
36 (
37  const primitiveMesh& mesh,
38  const pointField& p,
39  const vectorField& fCtrs,
40  const vectorField& fAreas,
41 
42  const label facei,
43  const point& ownCc,
44  const point& neiCc
45 )
46 {
47  vector Cpf = fCtrs[facei] - ownCc;
48  vector d = neiCc - ownCc;
49 
50  // Skewness vector
51  vector sv =
52  Cpf
53  - ((fAreas[facei] & Cpf)/((fAreas[facei] & d) + ROOTVSMALL))*d;
54  vector svHat = sv/(mag(sv) + ROOTVSMALL);
55 
56  // Normalisation distance calculated as the approximate distance
57  // from the face centre to the edge of the face in the direction
58  // of the skewness
59  scalar fd = 0.2*mag(d) + ROOTVSMALL;
60  const face& f = mesh.faces()[facei];
61  forAll(f, pi)
62  {
63  fd = max(fd, mag(svHat & (p[f[pi]] - fCtrs[facei])));
64  }
65 
66  // Normalised skewness
67  return mag(sv)/fd;
68 }
69 
70 
72 (
73  const primitiveMesh& mesh,
74  const pointField& p,
75  const vectorField& fCtrs,
76  const vectorField& fAreas,
77 
78  const label facei,
79  const point& ownCc
80 )
81 {
82  vector Cpf = fCtrs[facei] - ownCc;
83 
84  vector normal = normalised(fAreas[facei]);
85  vector d = normal*(normal & Cpf);
86 
87 
88  // Skewness vector
89  vector sv =
90  Cpf
91  - ((fAreas[facei] & Cpf)/((fAreas[facei] & d) + ROOTVSMALL))*d;
92  vector svHat = sv/(mag(sv) + ROOTVSMALL);
93 
94  // Normalisation distance calculated as the approximate distance
95  // from the face centre to the edge of the face in the direction
96  // of the skewness
97  scalar fd = 0.4*mag(d) + ROOTVSMALL;
98  const face& f = mesh.faces()[facei];
99  forAll(f, pi)
100  {
101  fd = max(fd, mag(svHat & (p[f[pi]] - fCtrs[facei])));
102  }
103 
104  // Normalised skewness
105  return mag(sv)/fd;
106 }
107 
108 
110 (
111  const point& ownCc,
112  const point& neiCc,
113  const vector& s
114 )
115 {
116  vector d = neiCc - ownCc;
117 
118  return (d & s)/(mag(d)*mag(s) + ROOTVSMALL);
119 }
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
125 (
126  const primitiveMesh& mesh,
127  const vectorField& areas,
128  const vectorField& cc
129 )
130 {
131  const labelList& own = mesh.faceOwner();
132  const labelList& nei = mesh.faceNeighbour();
133 
134  tmp<scalarField> tortho(new scalarField(mesh.nInternalFaces()));
135  scalarField& ortho = tortho.ref();
136 
137  // Internal faces
138  forAll(nei, facei)
139  {
140  ortho[facei] = faceOrthogonality
141  (
142  cc[own[facei]],
143  cc[nei[facei]],
144  areas[facei]
145  );
146  }
147 
148  return tortho;
149 }
150 
151 
153 (
154  const primitiveMesh& mesh,
155  const pointField& p,
156  const vectorField& fCtrs,
157  const vectorField& fAreas,
158  const vectorField& cellCtrs
159 )
160 {
161  const labelList& own = mesh.faceOwner();
162  const labelList& nei = mesh.faceNeighbour();
163 
164  tmp<scalarField> tskew(new scalarField(mesh.nFaces()));
165  scalarField& skew = tskew.ref();
166 
167  forAll(nei, facei)
168  {
169  skew[facei] = faceSkewness
170  (
171  mesh,
172  p,
173  fCtrs,
174  fAreas,
175 
176  facei,
177  cellCtrs[own[facei]],
178  cellCtrs[nei[facei]]
179  );
180  }
181 
182 
183  // Boundary faces: consider them to have only skewness error.
184  // (i.e. treat as if mirror cell on other side)
185 
186  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
187  {
188  skew[facei] = boundaryFaceSkewness
189  (
190  mesh,
191  p,
192  fCtrs,
193  fAreas,
194  facei,
195  cellCtrs[own[facei]]
196  );
197  }
198 
199  return tskew;
200 }
201 
202 
204 (
205  const primitiveMesh& mesh,
206  const pointField& points,
207  const vectorField& ctrs,
208 
209  scalarField& ownPyrVol,
210  scalarField& neiPyrVol
211 )
212 {
213  const labelList& own = mesh.faceOwner();
214  const labelList& nei = mesh.faceNeighbour();
215  const faceList& f = mesh.faces();
216 
217  ownPyrVol.setSize(mesh.nFaces());
218  neiPyrVol.setSize(mesh.nInternalFaces());
219 
220  forAll(f, facei)
221  {
222  // Create the owner pyramid
223  ownPyrVol[facei] = -pyramidPointFaceRef
224  (
225  f[facei],
226  ctrs[own[facei]]
227  ).mag(points);
228 
229  if (mesh.isInternalFace(facei))
230  {
231  // Create the neighbour pyramid - it will have positive volume
232  neiPyrVol[facei] = pyramidPointFaceRef
233  (
234  f[facei],
235  ctrs[nei[facei]]
236  ).mag(points);
237  }
238  }
239 }
240 
241 
243 (
244  const primitiveMesh& mesh,
245  const Vector<label>& meshD,
246  const vectorField& areas,
247  const scalarField& vols,
248 
249  scalarField& openness,
250  scalarField& aratio
251 )
252 {
253  const labelList& own = mesh.faceOwner();
254  const labelList& nei = mesh.faceNeighbour();
255 
256  // Loop through cell faces and sum up the face area vectors for each cell.
257  // This should be zero in all vector components
258 
259  vectorField sumClosed(mesh.nCells(), Zero);
260  vectorField sumMagClosed(mesh.nCells(), Zero);
261 
262  forAll(own, facei)
263  {
264  // Add to owner
265  sumClosed[own[facei]] += areas[facei];
266  sumMagClosed[own[facei]] += cmptMag(areas[facei]);
267  }
268 
269  forAll(nei, facei)
270  {
271  // Subtract from neighbour
272  sumClosed[nei[facei]] -= areas[facei];
273  sumMagClosed[nei[facei]] += cmptMag(areas[facei]);
274  }
275 
276 
277  label nDims = 0;
278  for (direction dir = 0; dir < vector::nComponents; dir++)
279  {
280  if (meshD[dir] == 1)
281  {
282  nDims++;
283  }
284  }
285 
286 
287  // Check the sums
288  openness.setSize(mesh.nCells());
289  aratio.setSize(mesh.nCells());
290 
291  forAll(sumClosed, celli)
292  {
293  scalar maxOpenness = 0;
294 
295  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
296  {
297  maxOpenness = max
298  (
299  maxOpenness,
300  mag(sumClosed[celli][cmpt])
301  /(sumMagClosed[celli][cmpt] + ROOTVSMALL)
302  );
303  }
304  openness[celli] = maxOpenness;
305 
306  // Calculate the aspect ration as the maximum of Cartesian component
307  // aspect ratio to the total area hydraulic area aspect ratio
308  scalar minCmpt = VGREAT;
309  scalar maxCmpt = -VGREAT;
310  for (direction dir = 0; dir < vector::nComponents; dir++)
311  {
312  if (meshD[dir] == 1)
313  {
314  minCmpt = min(minCmpt, sumMagClosed[celli][dir]);
315  maxCmpt = max(maxCmpt, sumMagClosed[celli][dir]);
316  }
317  }
318 
319  scalar aspectRatio = maxCmpt/(minCmpt + ROOTVSMALL);
320  if (nDims == 3)
321  {
322  scalar v = max(ROOTVSMALL, vols[celli]);
323 
324  aspectRatio = max
325  (
326  aspectRatio,
327  1.0/6.0*cmptSum(sumMagClosed[celli])/pow(v, 2.0/3.0)
328  );
329  }
330 
331  aratio[celli] = aspectRatio;
332  }
333 }
334 
335 
337 (
338  const scalar maxSin,
339  const primitiveMesh& mesh,
340  const pointField& p,
341  const vectorField& faceAreas
342 )
343 {
344  const faceList& fcs = mesh.faces();
345 
346  vectorField faceNormals(faceAreas);
347  faceNormals /= mag(faceNormals) + ROOTVSMALL;
348 
349  tmp<scalarField> tfaceAngles(new scalarField(mesh.nFaces()));
350  scalarField& faceAngles = tfaceAngles.ref();
351 
352 
353  forAll(fcs, facei)
354  {
355  const face& f = fcs[facei];
356 
357  // Get edge from f[0] to f[size-1];
358  vector ePrev(p[f.first()] - p[f.last()]);
359  scalar magEPrev = mag(ePrev);
360  ePrev /= magEPrev + ROOTVSMALL;
361 
362  scalar maxEdgeSin = 0.0;
363 
364  forAll(f, fp0)
365  {
366  // Get vertex after fp
367  label fp1 = f.fcIndex(fp0);
368 
369  // Normalized vector between two consecutive points
370  vector e10(p[f[fp1]] - p[f[fp0]]);
371  scalar magE10 = mag(e10);
372  e10 /= magE10 + ROOTVSMALL;
373 
374  if (magEPrev > SMALL && magE10 > SMALL)
375  {
376  vector edgeNormal = ePrev ^ e10;
377  scalar magEdgeNormal = mag(edgeNormal);
378 
379  if (magEdgeNormal < maxSin)
380  {
381  // Edges (almost) aligned -> face is ok.
382  }
383  else
384  {
385  // Check normal
386  edgeNormal /= magEdgeNormal;
387 
388  if ((edgeNormal & faceNormals[facei]) < SMALL)
389  {
390  maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
391  }
392  }
393  }
394 
395  ePrev = e10;
396  magEPrev = magE10;
397  }
398 
399  faceAngles[facei] = maxEdgeSin;
400  }
401 
402  return tfaceAngles;
403 }
404 
405 
407 (
408  const primitiveMesh& mesh,
409  const pointField& p,
410  const vectorField& fCtrs,
411  const vectorField& faceAreas
412 )
413 {
414  const faceList& fcs = mesh.faces();
415 
416  // Areas are calculated as the sum of areas. (see
417  // primitiveMeshFaceCentresAndAreas.C)
418  scalarField magAreas(mag(faceAreas));
419 
420  tmp<scalarField> tfaceFlatness(new scalarField(mesh.nFaces(), 1.0));
421  scalarField& faceFlatness = tfaceFlatness.ref();
422 
423  typedef Vector<solveScalar> solveVector;
424 
425  forAll(fcs, facei)
426  {
427  const face& f = fcs[facei];
428 
429  if (f.size() > 3 && magAreas[facei] > ROOTVSMALL)
430  {
431  const solveVector fc = fCtrs[facei];
432 
433  // Calculate the sum of magnitude of areas and compare to magnitude
434  // of sum of areas.
435 
436  solveScalar sumA = 0.0;
437 
438  forAll(f, fp)
439  {
440  const solveVector thisPoint = p[f[fp]];
441  const solveVector nextPoint = p[f.nextLabel(fp)];
442 
443  // Triangle around fc.
444  solveVector n = 0.5*((nextPoint - thisPoint)^(fc - thisPoint));
445  sumA += mag(n);
446  }
447 
448  faceFlatness[facei] = magAreas[facei]/(sumA + ROOTVSMALL);
449  }
450  }
451 
452  return tfaceFlatness;
453 }
454 
455 
457 (
458  const primitiveMesh& mesh,
459  const Vector<label>& meshD,
460  const vectorField& faceAreas,
461  const bitSet& internalOrCoupledFace
462 )
463 {
464  // Determine number of dimensions and (for 2D) missing dimension
465  label nDims = 0;
466  label twoD = -1;
467  for (direction dir = 0; dir < vector::nComponents; dir++)
468  {
469  if (meshD[dir] == 1)
470  {
471  nDims++;
472  }
473  else
474  {
475  twoD = dir;
476  }
477  }
478 
479  tmp<scalarField> tcellDeterminant(new scalarField(mesh.nCells()));
480  scalarField& cellDeterminant = tcellDeterminant.ref();
481 
482  const cellList& c = mesh.cells();
483 
484  if (nDims == 1)
485  {
486  cellDeterminant = 1.0;
487  }
488  else
489  {
490  forAll(c, celli)
491  {
492  const labelList& curFaces = c[celli];
493 
494  // Calculate local normalization factor
495  scalar avgArea = 0;
496 
497  label nInternalFaces = 0;
498 
499  forAll(curFaces, i)
500  {
501  if (internalOrCoupledFace.test(curFaces[i]))
502  {
503  avgArea += mag(faceAreas[curFaces[i]]);
504 
505  nInternalFaces++;
506  }
507  }
508 
509  if (nInternalFaces == 0 || avgArea < ROOTVSMALL)
510  {
511  cellDeterminant[celli] = 0;
512  }
513  else
514  {
515  avgArea /= nInternalFaces;
516 
517  symmTensor areaTensor(Zero);
518 
519  forAll(curFaces, i)
520  {
521  if (internalOrCoupledFace.test(curFaces[i]))
522  {
523  areaTensor += sqr(faceAreas[curFaces[i]]/avgArea);
524  }
525  }
526 
527  if (nDims == 2)
528  {
529  // Add the missing eigenvector (such that it does not
530  // affect the determinant)
531  if (twoD == 0)
532  {
533  areaTensor.xx() = 1;
534  }
535  else if (twoD == 1)
536  {
537  areaTensor.yy() = 1;
538  }
539  else
540  {
541  areaTensor.zz() = 1;
542  }
543  }
544 
545  // Note:
546  // - normalise to be 0..1 (since cube has eigenvalues 2 2 2)
547  // - we use the determinant (i.e. 3rd invariant) and not e.g.
548  // condition number (= max ev / min ev) since we are
549  // interested in the minimum connectivity and not the
550  // uniformity. Using the condition number on corner cells
551  // leads to uniformity 1 i.e. equally bad in all three
552  // directions which is not what we want.
553  cellDeterminant[celli] = mag(det(areaTensor))/8.0;
554  }
555  }
556  }
557 
558  return tcellDeterminant;
559 }
560 
561 
562 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::primitiveMeshTools::cellDeterminant
static tmp< scalarField > cellDeterminant(const primitiveMesh &mesh, const Vector< label > &directions, const vectorField &faceAreas, const bitSet &internalOrCoupledFace)
Generate cell determinant field. Normalised to 1 for an internal cube.
Definition: primitiveMeshTools.C:457
Foam::SymmTensor< scalar >
p
volScalarField & p
Definition: createFieldRefs.H:8
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:64
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:138
primitiveMeshTools.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::primitiveMeshTools::faceConcavity
static tmp< scalarField > faceConcavity(const scalar maxSin, const primitiveMesh &mesh, const pointField &p, const vectorField &faceAreas)
Generate face concavity field. Returns per face the (sin of the)
Definition: primitiveMeshTools.C:337
Foam::primitiveMeshTools::faceOrthogonality
static tmp< scalarField > faceOrthogonality(const primitiveMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate non-orthogonality field (internal faces only)
Definition: primitiveMeshTools.C:125
Foam::bitSet::test
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:512
syncTools.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
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::SymmTensor::xx
const Cmpt & xx() const
Definition: SymmTensorI.H:84
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
faceNormals
surfaceVectorField faceNormals(mesh.Sf()/mesh.magSf())
Foam::Field< vector >
Foam::primitiveMeshTools::faceFlatness
static tmp< scalarField > faceFlatness(const primitiveMesh &mesh, const pointField &p, const vectorField &fCtrs, const vectorField &faceAreas)
Generate face flatness field. Compares the individual triangles'.
Definition: primitiveMeshTools.C:407
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::cmptSum
Cmpt cmptSum(const SphericalTensor< Cmpt > &st)
Return the sum of components of a SphericalTensor.
Definition: SphericalTensorI.H:164
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:483
Foam::primitiveMeshTools::faceSkewness
static tmp< scalarField > faceSkewness(const primitiveMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
Definition: primitiveMeshTools.C:153
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::SymmTensor::yy
const Cmpt & yy() const
Definition: SymmTensorI.H:108
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
Foam::primitiveMeshTools::cellClosedness
static void cellClosedness(const primitiveMesh &mesh, const Vector< label > &meshD, const vectorField &areas, const scalarField &vols, scalarField &openness, scalarField &aratio)
Generate cell openness and cell ascpect ratio field.
Definition: primitiveMeshTools.C:243
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::primitiveMeshTools::boundaryFaceSkewness
static scalar boundaryFaceSkewness(const primitiveMesh &mesh, const pointField &p, const vectorField &fCtrs, const vectorField &fAreas, const label facei, const point &ownCc)
Skewness of single boundary face.
Definition: primitiveMeshTools.C:72
Foam::direction
uint8_t direction
Definition: direction.H:47
Foam::SymmTensor::zz
const Cmpt & zz() const
Definition: SymmTensorI.H:132
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::det
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:62
Foam::pyramidPointFaceRef
pyramid< point, const point &, const face & > pyramidPointFaceRef
Definition: pyramidPointFaceRef.H:47
pyramidPointFaceRef.H
Foam::primitiveMeshTools::facePyramidVolume
static void facePyramidVolume(const primitiveMesh &mesh, const pointField &points, const vectorField &cellCtrs, scalarField &ownPyrVol, scalarField &neiPyrVol)
Generate face pyramid volume fields.
Definition: primitiveMeshTools.C:204
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78