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