Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "primitiveMesh.H"
31#include "syncTools.H"
32#include "pyramid.H"
33#include "tetrahedron.H"
34#include "PrecisionAdaptor.H"
35
36// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41// Average of points
42// Note: use double precision to avoid overflows when summing
44(
45 const UList<point>& points,
47)
48{
49 doubleVector avg(Zero);
50
51 if (const auto n = pointLabels.size(); n)
52 {
53 for (const auto pointi : pointLabels)
54 {
55 avg += points[pointi];
56 }
57 avg /= n;
58 }
59
60 return avg;
62
63} // End namespace Foam
64
65
66// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
67
69(
70 const primitiveMesh& mesh,
71 const UList<label>& faceIDs,
72 const pointField& p,
73 vectorField& fCtrs,
74 vectorField& fAreas
75)
76{
77 const faceList& fs = mesh.faces();
78
79 for (const label facei : faceIDs)
80 {
81 const auto& f = fs[facei];
82 const label nPoints = f.size();
83
84 // If the face is a triangle, do a direct calculation for efficiency
85 // and to avoid round-off error-related problems
86 if (nPoints == 3)
87 {
88 const triPointRef tri(p, f[0], f[1], f[2]);
89 fCtrs[facei] = tri.centre();
90 fAreas[facei] = tri.areaNormal();
91 }
92 else
93 {
94 solveVector sumN(Zero);
95 solveScalar sumA(0);
96 solveVector sumAc(Zero);
97
98 // Estimated centre by averaging the face points
99 const solveVector fCentre(pointsAverage(p, f));
100
101 for (label pi = 0; pi < nPoints; ++pi)
102 {
103 const label nextPi(pi == nPoints-1 ? 0 : pi+1);
104 const solveVector nextPoint(p[f[nextPi]]);
105 const solveVector thisPoint(p[f[pi]]);
106
107 solveVector c = thisPoint + nextPoint + fCentre;
108 solveVector n = (nextPoint - thisPoint)^(fCentre - thisPoint);
109 solveScalar a = mag(n);
110 sumN += n;
111 sumA += a;
112 sumAc += a*c;
113 }
114
115 // This is to deal with zero-area faces. Mark very small faces
116 // to be detected in e.g., processorPolyPatch.
117 if (sumA < ROOTVSMALL)
118 {
119 fCtrs[facei] = fCentre;
120 fAreas[facei] = Zero;
121 }
122 else
123 {
124 fCtrs[facei] = (1.0/3.0)*sumAc/sumA;
125 fAreas[facei] = 0.5*sumN;
126 }
127 }
128 }
129}
130
131
133(
134 const primitiveMesh& mesh,
135 const vectorField& fCtrs,
136 const vectorField& fAreas,
137 const List<label>& cellIDs,
138 const List<cell>& cells,
139 vectorField& cellCtrs_s,
140 scalarField& cellVols_s
141)
142{
143 PrecisionAdaptor<solveVector, vector> tcellCtrs(cellCtrs_s, false);
144 PrecisionAdaptor<solveScalar, scalar> tcellVols(cellVols_s, false);
145 auto& cellCtrs = tcellCtrs.ref();
146 auto& cellVols = tcellVols.ref();
147
148
149 // Use current cell centres as estimates for the new cell centres
150 // Note - not currently used
151 // - Introduces a small difference (seen when write precision extended to
152 // 16) to the cell centre and volume values, typically last 4 digits
153 //const vectorField Cc0(cellCtrs, cellIDs);
154
155 // Estimate cell centres using current face centres
156 // - Same method as used by makeCellCentresAndVols()
157 vectorField Cc0(cellIDs.size(), Zero);
158 {
159 labelField nCellFaces(cellIDs.size(), Zero);
160 const auto& cells = mesh.cells();
161
162 forAll(cellIDs, i)
163 {
164 const label celli = cellIDs[i];
165 const cell& c = cells[celli];
166 for (const auto facei : c)
167 {
168 Cc0[i] += fCtrs[facei];
169 ++nCellFaces[i];
170 }
171 }
172
173 forAll(Cc0, i)
174 {
175 Cc0[i] /= nCellFaces[i];
176 }
177 }
178
179 // Clear the fields for accumulation
180 for (const label celli : cellIDs)
181 {
182 cellCtrs[celli] = Zero;
183 cellVols[celli] = Zero;
184 }
185
186 const auto& own = mesh.faceOwner();
187
188 forAll(cellIDs, i)
189 {
190 const label celli = cellIDs[i];
191 const auto& c = cells[celli];
192 const solveVector cc(Cc0[i]);
193
194 for (const label facei : c)
195 {
196 const solveVector fc(fCtrs[facei]);
197 const solveVector fA(fAreas[facei]);
198
199 const solveScalar pyr3Vol = own[facei] == celli
200 ? fA & (fc - cc)
201 : fA & (cc - fc);
202
203 // Calculate face-pyramid centre
204 const solveVector pc = (3.0/4.0)*fc + (1.0/4.0)*cc;
205
206 // Accumulate volume-weighted face-pyramid centre
207 cellCtrs[celli] += pyr3Vol*pc;
208
209 // Accumulate face-pyramid volume
210 cellVols[celli] += pyr3Vol;
211 }
212
213 if (mag(cellVols[celli]) > VSMALL)
214 {
215 cellCtrs[celli] /= cellVols[celli];
216 cellVols[celli] *= (1.0/3.0);
217 }
218 else
220 cellCtrs[celli] = Cc0[i];
221 }
222 }
223}
224
225
227(
228 const UList<face>& fcs,
229 const pointField& p,
230 vectorField& fCtrs,
231 vectorField& fAreas
232)
233{
234 // Safety first - ensure properly sized
235 fCtrs.resize_nocopy(fcs.size());
236 fAreas.resize_nocopy(fcs.size());
237
238 forAll(fcs, facei)
239 {
240 const face& f = fcs[facei];
241 const label nPoints = f.size();
242
243 // If the face is a triangle, do a direct calculation for efficiency
244 // and to avoid round-off error-related problems
245 if (nPoints == 3)
246 {
247 const triPointRef tri(p, f[0], f[1], f[2]);
248 fCtrs[facei] = tri.centre();
249 fAreas[facei] = tri.areaNormal();
250 }
251 else
252 {
253 solveVector sumN(Zero);
254 solveScalar sumA(0);
255 solveVector sumAc(Zero);
256
257 // Estimated centre by averaging the face points
258 const solveVector fCentre(pointsAverage(p, f));
259
260 for (label pi = 0; pi < nPoints; ++pi)
261 {
262 const solveVector thisPoint(p[f.thisLabel(pi)]);
263 const solveVector nextPoint(p[f.nextLabel(pi)]);
264
265 solveVector c = thisPoint + nextPoint + fCentre;
266 solveVector n = (nextPoint - thisPoint)^(fCentre - thisPoint);
267 solveScalar a = mag(n);
268
269 sumN += n;
270 sumA += a;
271 sumAc += a*c;
272 }
273
274 // This is to deal with zero-area faces. Mark very small faces
275 // to be detected in e.g., processorPolyPatch.
276 if (sumA < ROOTVSMALL)
277 {
278 fCtrs[facei] = fCentre;
279 fAreas[facei] = Zero;
280 }
281 else
282 {
283 fCtrs[facei] = (1.0/3.0)*sumAc/sumA;
284 fAreas[facei] = 0.5*sumN;
285 }
286 }
287 }
288}
289
290
292(
293 const primitiveMesh& mesh,
294 const pointField& p,
295 vectorField& fCtrs,
297)
298{
299 makeFaceCentresAndAreas(mesh.faces(), p, fCtrs, fAreas);
300}
301
302
304(
305 const primitiveMesh& mesh,
306 const vectorField& fCtrs,
307 const vectorField& fAreas,
308 vectorField& cellCtrs_s,
309 scalarField& cellVols_s
310)
311{
312 PrecisionAdaptor<solveVector, vector> tcellCtrs(cellCtrs_s, false);
313 PrecisionAdaptor<solveScalar, scalar> tcellVols(cellVols_s, false);
314 auto& cellCtrs = tcellCtrs.ref();
315 auto& cellVols = tcellVols.ref();
316
317 // Clear the fields for accumulation
318 cellCtrs = Zero;
319 cellVols = Zero;
320
321 const labelList& own = mesh.faceOwner();
322 const labelList& nei = mesh.faceNeighbour();
323
324 // first estimate the approximate cell centre as the average of
325 // face centres
326
327 Field<solveVector> cEst(mesh.nCells(), Zero);
328 labelField nCellFaces(mesh.nCells(), Zero);
329
330 forAll(own, facei)
331 {
332 cEst[own[facei]] += fCtrs[facei];
333 ++nCellFaces[own[facei]];
334 }
335
336 forAll(nei, facei)
337 {
338 cEst[nei[facei]] += fCtrs[facei];
339 ++nCellFaces[nei[facei]];
340 }
341
342 forAll(cEst, celli)
343 {
344 cEst[celli] /= nCellFaces[celli];
345 }
346
347 forAll(own, facei)
348 {
349 const solveVector fc(fCtrs[facei]);
350 const solveVector fA(fAreas[facei]);
351
352 // Calculate 3*face-pyramid volume
353 solveScalar pyr3Vol = fA & (fc - cEst[own[facei]]);
354
355 // Calculate face-pyramid centre
356 solveVector pc = (3.0/4.0)*fc + (1.0/4.0)*cEst[own[facei]];
357
358 // Accumulate volume-weighted face-pyramid centre
359 cellCtrs[own[facei]] += pyr3Vol*pc;
360
361 // Accumulate face-pyramid volume
362 cellVols[own[facei]] += pyr3Vol;
363 }
364
365 forAll(nei, facei)
366 {
367 const solveVector fc(fCtrs[facei]);
368 const solveVector fA(fAreas[facei]);
369
370 // Calculate 3*face-pyramid volume
371 solveScalar pyr3Vol = fA & (cEst[nei[facei]] - fc);
372
373 // Calculate face-pyramid centre
374 solveVector pc = (3.0/4.0)*fc + (1.0/4.0)*cEst[nei[facei]];
375
376 // Accumulate volume-weighted face-pyramid centre
377 cellCtrs[nei[facei]] += pyr3Vol*pc;
378
379 // Accumulate face-pyramid volume
380 cellVols[nei[facei]] += pyr3Vol;
381 }
382
383 forAll(cellCtrs, celli)
384 {
385 if (mag(cellVols[celli]) > VSMALL)
386 {
387 cellCtrs[celli] /= cellVols[celli];
388 }
389 else
390 {
391 cellCtrs[celli] = cEst[celli];
393 }
394
395 cellVols *= (1.0/3.0);
396}
397
398
400(
401 const UList<face>& fcs,
402 const pointField& p,
403 const vectorField& fCtrs,
404 const vectorField& fAreas,
405
406 const label facei,
407 const point& ownCc,
408 const point& neiCc
409)
410{
411 const face& f = fcs[facei];
412
413 vector Cpf = fCtrs[facei] - ownCc;
414 vector d = neiCc - ownCc;
415
416 // Skewness vector
417 vector sv =
418 Cpf
419 - ((fAreas[facei] & Cpf)/((fAreas[facei] & d) + ROOTVSMALL))*d;
420 vector svHat = sv/(mag(sv) + ROOTVSMALL);
421
422 // Normalisation distance calculated as the approximate distance
423 // from the face centre to the edge of the face in the direction
424 // of the skewness
425 scalar fd = 0.2*mag(d) + ROOTVSMALL;
426 forAll(f, pi)
427 {
428 fd = max(fd, mag(svHat & (p[f[pi]] - fCtrs[facei])));
430
431 // Normalised skewness
432 return mag(sv)/fd;
433}
434
435
437(
438 const UList<face>& fcs,
439 const pointField& p,
440 const vectorField& fCtrs,
441 const vectorField& fAreas,
442
443 const label facei,
444 const point& ownCc
445)
446{
447 const face& f = fcs[facei];
448
449 vector Cpf = fCtrs[facei] - ownCc;
450
451 vector normal = normalised(fAreas[facei]);
452 vector d = normal*(normal & Cpf);
453
454 // Skewness vector
455 vector sv =
456 Cpf
457 - ((fAreas[facei] & Cpf)/((fAreas[facei] & d) + ROOTVSMALL))*d;
458 vector svHat = sv/(mag(sv) + ROOTVSMALL);
459
460 // Normalisation distance calculated as the approximate distance
461 // from the face centre to the edge of the face in the direction
462 // of the skewness
463 scalar fd = 0.4*mag(d) + ROOTVSMALL;
464 forAll(f, pi)
465 {
466 fd = max(fd, mag(svHat & (p[f[pi]] - fCtrs[facei])));
468
469 // Normalised skewness
470 return mag(sv)/fd;
471}
472
473
475(
476 const primitiveMesh& mesh,
477 const pointField& p,
478 const vectorField& fCtrs,
479 const vectorField& fAreas,
480
481 const label facei,
482 const point& ownCc,
483 const point& neiCc
484)
485{
486 return faceSkewness(mesh.faces(), p, fCtrs, fAreas, facei, ownCc, neiCc);
487}
488
489
491(
492 const primitiveMesh& mesh,
493 const pointField& p,
494 const vectorField& fCtrs,
495 const vectorField& fAreas,
496
497 const label facei,
498 const point& ownCc
499)
500{
501 return boundaryFaceSkewness(mesh.faces(), p, fCtrs, fAreas, facei, ownCc);
502}
503
504
506(
507 const point& ownCc,
508 const point& neiCc,
509 const vector& s
510)
511{
512 vector d = neiCc - ownCc;
514 return (d & s)/(mag(d)*mag(s) + ROOTVSMALL);
515}
516
517
518// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
519
521(
522 const primitiveMesh& mesh,
523 const vectorField& areas,
524 const vectorField& cc
525)
526{
527 const labelList& own = mesh.faceOwner();
528 const labelList& nei = mesh.faceNeighbour();
529
530 auto tortho = tmp<scalarField>::New(mesh.nInternalFaces());
531 auto& ortho = tortho.ref();
532
533 // Internal faces
534 forAll(nei, facei)
535 {
536 ortho[facei] = faceOrthogonality
537 (
538 cc[own[facei]],
539 cc[nei[facei]],
540 areas[facei]
541 );
542 }
543
544 return tortho;
545}
546
547
549(
550 const primitiveMesh& mesh,
551 const pointField& p,
552 const vectorField& fCtrs,
553 const vectorField& fAreas,
554 const vectorField& cellCtrs
555)
556{
557 const labelList& own = mesh.faceOwner();
558 const labelList& nei = mesh.faceNeighbour();
559 const faceList& faces = mesh.faces();
560
561 auto tskew = tmp<scalarField>::New(mesh.nFaces());
562 auto& skew = tskew.ref();
563
564 // Internal faces
565 for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
566 {
567 skew[facei] = faceSkewness
568 (
569 faces,
570 p,
571 fCtrs,
572 fAreas,
573
574 facei,
575 cellCtrs[own[facei]],
576 cellCtrs[nei[facei]]
577 );
578 }
579
580
581 // Boundary faces: consider them to have only skewness error.
582 // (i.e. treat as if mirror cell on other side)
583
584 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
585 {
586 skew[facei] = boundaryFaceSkewness
587 (
588 faces,
589 p,
590 fCtrs,
591 fAreas,
592 facei,
593 cellCtrs[own[facei]]
594 );
595 }
596
597 return tskew;
598}
599
600
602(
603 const primitiveMesh& mesh,
604 const pointField& points,
605 const vectorField& ctrs,
606
607 scalarField& ownPyrVol,
608 scalarField& neiPyrVol
609)
610{
611 const labelList& own = mesh.faceOwner();
612 const labelList& nei = mesh.faceNeighbour();
613 const faceList& f = mesh.faces();
614
615 ownPyrVol.setSize(mesh.nFaces());
616 neiPyrVol.setSize(mesh.nInternalFaces());
617
618 forAll(f, facei)
619 {
620 // Create the owner pyramid
621 ownPyrVol[facei] = -pyramidPointFaceRef
622 (
623 f[facei],
624 ctrs[own[facei]]
625 ).mag(points);
626
627 if (mesh.isInternalFace(facei))
628 {
629 // Create the neighbour pyramid - it will have positive volume
630 neiPyrVol[facei] = pyramidPointFaceRef
631 (
632 f[facei],
633 ctrs[nei[facei]]
634 ).mag(points);
635 }
636 }
637}
638
639
641(
642 const primitiveMesh& mesh,
643 const Vector<label>& meshD,
644 const vectorField& areas,
645 const scalarField& vols,
646
647 scalarField& openness,
648 scalarField& aratio
649)
650{
651 const labelList& own = mesh.faceOwner();
652 const labelList& nei = mesh.faceNeighbour();
653
654 // Loop through cell faces and sum up the face area vectors for each cell.
655 // This should be zero in all vector components
656
657 vectorField sumClosed(mesh.nCells(), Zero);
658 vectorField sumMagClosed(mesh.nCells(), Zero);
659
660 forAll(own, facei)
661 {
662 // Add to owner
663 sumClosed[own[facei]] += areas[facei];
664 sumMagClosed[own[facei]] += cmptMag(areas[facei]);
665 }
666
667 forAll(nei, facei)
668 {
669 // Subtract from neighbour
670 sumClosed[nei[facei]] -= areas[facei];
671 sumMagClosed[nei[facei]] += cmptMag(areas[facei]);
672 }
673
674
675 label nDims = 0;
676 for (direction dir = 0; dir < vector::nComponents; dir++)
677 {
678 if (meshD[dir] == 1)
679 {
680 nDims++;
681 }
682 }
683
684
685 // Check the sums
686 openness.setSize(mesh.nCells());
687 aratio.setSize(mesh.nCells());
688
689 forAll(sumClosed, celli)
690 {
691 scalar maxOpenness = 0;
692
693 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
694 {
695 maxOpenness = max
696 (
697 maxOpenness,
698 mag(sumClosed[celli][cmpt])
699 /(sumMagClosed[celli][cmpt] + ROOTVSMALL)
700 );
701 }
702 openness[celli] = maxOpenness;
703
704 // Calculate the aspect ration as the maximum of Cartesian component
705 // aspect ratio to the total area hydraulic area aspect ratio
706 scalar minCmpt = VGREAT;
707 scalar maxCmpt = -VGREAT;
708 for (direction dir = 0; dir < vector::nComponents; dir++)
709 {
710 if (meshD[dir] == 1)
711 {
712 minCmpt = min(minCmpt, sumMagClosed[celli][dir]);
713 maxCmpt = max(maxCmpt, sumMagClosed[celli][dir]);
714 }
715 }
716
717 scalar aspectRatio = maxCmpt/(minCmpt + ROOTVSMALL);
718 if (nDims == 3)
719 {
720 scalar v = max(ROOTVSMALL, vols[celli]);
721
722 aspectRatio = max
723 (
724 aspectRatio,
725 1.0/6.0*cmptSum(sumMagClosed[celli])/pow(v, 2.0/3.0)
726 );
728
729 aratio[celli] = aspectRatio;
730 }
731}
732
733
734Foam::tmp<Foam::scalarField> Foam::primitiveMeshTools::faceConcavity
735(
736 const scalar maxSin,
737 const primitiveMesh& mesh,
738 const pointField& p,
739 const vectorField& faceAreas
740)
741{
742 const faceList& fcs = mesh.faces();
743
744 vectorField faceNormals(faceAreas);
745 faceNormals /= mag(faceNormals) + ROOTVSMALL;
746
747 auto tfaceAngles = tmp<scalarField>::New(mesh.nFaces());
748 auto&& faceAngles = tfaceAngles.ref();
749
750 forAll(fcs, facei)
751 {
752 const face& f = fcs[facei];
753
754 // Normalized vector from f[size-1] to f[0];
755 vector ePrev(p[f.first()] - p[f.last()]);
756 scalar magEPrev = mag(ePrev);
757 ePrev /= magEPrev + ROOTVSMALL;
758
759 scalar maxEdgeSin = 0.0;
760
761 forAll(f, fp0)
762 {
763 // Normalized vector between two consecutive points
764 vector e10(p[f.nextLabel(fp0)] - p[f.thisLabel(fp0)]);
765 scalar magE10 = mag(e10);
766 e10 /= magE10 + ROOTVSMALL;
767
768 if (magEPrev > SMALL && magE10 > SMALL)
769 {
770 vector edgeNormal = ePrev ^ e10;
771 scalar magEdgeNormal = mag(edgeNormal);
772
773 if (magEdgeNormal < maxSin)
774 {
775 // Edges (almost) aligned -> face is ok.
776 }
777 else
778 {
779 // Check normal
780 edgeNormal /= magEdgeNormal;
781
782 if ((edgeNormal & faceNormals[facei]) < SMALL)
783 {
784 maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
785 }
786 }
787 }
788
789 ePrev = e10;
790 magEPrev = magE10;
791 }
792
793 faceAngles[facei] = maxEdgeSin;
794 }
795
796 return tfaceAngles;
797}
798
799
801(
802 const primitiveMesh& mesh,
803 const pointField& p,
804 const vectorField& fCtrs,
805 const vectorField& faceAreas
806)
807{
808 const faceList& fcs = mesh.faces();
809
810 // Areas are calculated as the sum of areas. (see
811 // primitiveMeshFaceCentresAndAreas.C)
812 scalarField magAreas(mag(faceAreas));
813
814 auto tfaceFlatness = tmp<scalarField>::New(mesh.nFaces(), scalar(1));
815 auto& faceFlatness = tfaceFlatness.ref();
816
817 forAll(fcs, facei)
818 {
819 const face& f = fcs[facei];
820
821 if (f.size() > 3 && magAreas[facei] > ROOTVSMALL)
822 {
823 const solveVector fc = fCtrs[facei];
824
825 // Calculate the sum of magnitude of areas and compare to magnitude
826 // of sum of areas.
827
828 solveScalar sumA = 0.0;
829
830 forAll(f, fp)
831 {
832 const solveVector thisPoint = p[f[fp]];
833 const solveVector nextPoint = p[f.nextLabel(fp)];
834
835 // Triangle around fc.
836 solveVector n = 0.5*((nextPoint - thisPoint)^(fc - thisPoint));
837 sumA += mag(n);
838 }
839
840 faceFlatness[facei] = magAreas[facei]/(sumA + ROOTVSMALL);
842 }
843
844 return tfaceFlatness;
845}
846
847
848Foam::tmp<Foam::scalarField> Foam::primitiveMeshTools::cellDeterminant
849(
850 const primitiveMesh& mesh,
851 const Vector<label>& meshD,
852 const vectorField& faceAreas,
853 const bitSet& internalOrCoupledFace
854)
855{
856 // Determine number of dimensions and (for 2D) missing dimension
857 label nDims = 0;
858 label twoD = -1;
859 for (direction dir = 0; dir < vector::nComponents; dir++)
860 {
861 if (meshD[dir] == 1)
862 {
863 nDims++;
864 }
865 else
866 {
867 twoD = dir;
868 }
869 }
870
871 auto tcellDeterminant = tmp<scalarField>::New(mesh.nCells());
872 auto& cellDeterminant = tcellDeterminant.ref();
873
874 const cellList& c = mesh.cells();
875
876 if (nDims == 1)
877 {
878 cellDeterminant = 1.0;
879 }
880 else
881 {
882 forAll(c, celli)
883 {
884 const labelList& curFaces = c[celli];
885
886 // Calculate local normalization factor
887 scalar avgArea = 0;
888
889 label nInternalFaces = 0;
890
891 forAll(curFaces, i)
892 {
893 if (internalOrCoupledFace.test(curFaces[i]))
894 {
895 avgArea += mag(faceAreas[curFaces[i]]);
896
897 nInternalFaces++;
898 }
899 }
900
901 if (nInternalFaces == 0 || avgArea < ROOTVSMALL)
902 {
903 cellDeterminant[celli] = 0;
904 }
905 else
906 {
907 avgArea /= nInternalFaces;
908
909 symmTensor areaTensor(Zero);
910
911 forAll(curFaces, i)
912 {
913 if (internalOrCoupledFace.test(curFaces[i]))
914 {
915 areaTensor += sqr(faceAreas[curFaces[i]]/avgArea);
916 }
917 }
918
919 if (nDims == 2)
920 {
921 // Add the missing eigenvector (such that it does not
922 // affect the determinant)
923 if (twoD == 0)
924 {
925 areaTensor.xx() = 1;
926 }
927 else if (twoD == 1)
928 {
929 areaTensor.yy() = 1;
930 }
931 else
932 {
933 areaTensor.zz() = 1;
934 }
935 }
936
937 // Note:
938 // - normalise to be 0..1 (since cube has eigenvalues 2 2 2)
939 // - we use the determinant (i.e. 3rd invariant) and not e.g.
940 // condition number (= max ev / min ev) since we are
941 // interested in the minimum connectivity and not the
942 // uniformity. Using the condition number on corner cells
943 // leads to uniformity 1 i.e. equally bad in all three
944 // directions which is not what we want.
945 cellDeterminant[celli] = mag(det(areaTensor))/8.0;
946 }
947 }
948 }
949
950 return tcellDeterminant;
951}
952
953
954// ************************************************************************* //
constexpr scalar pi(M_PI)
label n
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
void setSize(label n)
Alias for resize().
Definition List.H:536
A non-const Field/List wrapper with possible data conversion.
const Cmpt & yy() const noexcept
Definition SymmTensor.H:154
const Cmpt & xx() const noexcept
Definition SymmTensor.H:150
const Cmpt & zz() const noexcept
Definition SymmTensor.H:158
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
T & first()
Access first element of the list, position [0].
Definition UList.H:957
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
T & last()
Access last element of the list, position [size()-1].
Definition UList.H:971
static constexpr direction nComponents
Number of components in this vector space.
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
static void makeFaceCentresAndAreas(const UList< face > &faces, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Calculate face centres and areas for specified faces.
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) most concave angle between two conse...
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 aspect ratio field.
static scalar boundaryFaceSkewness(const UList< face > &faces, const pointField &p, const vectorField &fCtrs, const vectorField &fAreas, const label facei, const point &ownCc)
Skewness of single boundary face.
static tmp< scalarField > faceSkewness(const primitiveMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
static void updateCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, const List< label > &cellIDs, const List< cell > &cells, vectorField &cellCtrs_s, scalarField &cellVols_s)
Update cell centres and volumes for the cells in the set cellIDs.
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.
static tmp< scalarField > faceOrthogonality(const primitiveMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate non-orthogonality field (internal faces only).
static void updateFaceCentresAndAreas(const primitiveMesh &mesh, const UList< label > &faceIDs, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Update face centres and areas for the faces in the set faceIDs.
static void facePyramidVolume(const primitiveMesh &mesh, const pointField &points, const vectorField &cellCtrs, scalarField &ownPyrVol, scalarField &neiPyrVol)
Generate face pyramid volume fields.
static tmp< scalarField > faceFlatness(const primitiveMesh &mesh, const pointField &p, const vectorField &fCtrs, const vectorField &faceAreas)
Generate face flatness field. Compares the individual triangles' normals against the face average nor...
static void makeCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols)
Calculate cell centres and volumes from face properties.
Cell-face mesh analysis engine.
label nInternalFaces() const noexcept
Number of internal faces.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const cellList & cells() const
scalar mag(const UList< point > &points) const
Return scalar magnitude - returns volume of pyramid.
Definition pyramidI.H:69
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
static vector areaNormal(const Point &p0, const Point &p1, const Point &p2)
The area normal for a triangle defined by three points (right-hand rule). Magnitude equal to area of ...
Definition triangleI.H:169
static Point centre(const Point &p0, const Point &p1, const Point &p2)
The centre (centroid) of three points.
Definition triangleI.H:144
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
dynamicFvMesh & mesh
const pointField & points
label nPoints
const cellShapeList & cells
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))
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
static Foam::doubleVector pointsAverage(const UList< point > &points, const labelUList &pointLabels)
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Cmpt cmptSum(const SphericalTensor< Cmpt > &st)
Return the sum of components of a SphericalTensor.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Vector< double > doubleVector
Definition vector.H:54
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
pyramid< point, const point &, const face & > pyramidPointFaceRef
A pyramid using referred point and face.
Definition pyramid.H:70
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
Vector< solveScalar > solveVector
Definition vector.H:60
dimensionedTensor skew(const dimensionedTensor &dt)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const scalarField & cellVols