averageNeighbourFvGeometryScheme.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) 2020-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
30#include "fvMesh.H"
31#include "cellAspectRatio.H"
32#include "syncTools.H"
33#include "polyMeshTools.H"
34#include "unitConversion.H"
35#include "OBJstream.H"
36#include "surfaceWriter.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
44 (
47 dict
48 );
49}
50
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
55(
56 const scalar minRatio,
57 const vectorField& faceCentres,
59
60 vectorField& faceCorrection
61) const
62{
63 // Clip correction vector if any triangle becomes too small. Return number
64 // of correction vectors clipped
65
66 const pointField& p = mesh_.points();
67
68 label nClipped = 0;
69 for (label facei = 0; facei < mesh_.nFaces(); facei++)
70 {
71 #ifdef WM_SPDP
72 const solveVector fcCorr(faceCorrection[facei]);
73 #else
74 const vector& fcCorr = faceCorrection[facei];
75 #endif
76 if (fcCorr != solveVector::zero)
77 {
78 #ifdef WM_SPDP
79 const solveVector fn(faceNormals[facei]);
80 const solveVector fc(faceCentres[facei]);
81 #else
82 const vector& fn = faceNormals[facei];
83 const point& fc = faceCentres[facei];
84 #endif
85 const face& f = mesh_.faces()[facei];
86
87 forAll(f, fp)
88 {
89 const solveVector thisPt(p[f[fp]]);
90 const solveVector nextPt(p[f.fcValue(fp)]);
91 const solveVector d(nextPt-thisPt);
92
93 // Calculate triangle area with correction
94 const solveVector nCorr(d^(fc+fcCorr - thisPt));
95
96 if ((nCorr & fn) < 0)
97 {
98 // Triangle points wrong way
99 faceCorrection[facei] = vector::zero;
100 nClipped++;
101 break;
102 }
103 else
104 {
105 // Calculate triangle area without correction
106 const solveVector n(d^(fc - thisPt));
107 if ((n & fn) < 0)
108 {
109 // Original triangle points the wrong way, new one is ok
110 }
111 else
112 {
113 // Both point correctly. Make sure triangle doesn't get
114 // too small
115 if (mag(nCorr) < minRatio*mag(n))
116 {
117 faceCorrection[facei] = vector::zero;
118 nClipped++;
119 break;
120 }
121 }
122 }
123 }
124 }
125 }
126 return returnReduce(nClipped, sumOp<label>());
127}
128
129
131(
132 const pointField& cellCentres,
133 const vectorField& faceCentres,
135
136 scalarField& ownHeight,
137 scalarField& neiHeight
138) const
139{
140 ownHeight.setSize(mesh_.nFaces());
141 neiHeight.setSize(mesh_.nInternalFaces());
142
143 const labelList& own = mesh_.faceOwner();
144 const labelList& nei = mesh_.faceNeighbour();
145
146 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
147 {
148 const solveVector n = faceNormals[facei];
149 const solveVector fc = faceCentres[facei];
150 const solveVector ownCc = cellCentres[own[facei]];
151 const solveVector neiCc = cellCentres[nei[facei]];
152
153 ownHeight[facei] = ((fc-ownCc)&n);
154 neiHeight[facei] = ((neiCc-fc)&n);
155 }
156
157 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
158 {
159 const solveVector n = faceNormals[facei];
160 const solveVector fc = faceCentres[facei];
161 const solveVector ownCc = cellCentres[own[facei]];
162
163 ownHeight[facei] = ((fc-ownCc)&n);
164 }
165}
166
167
169(
170 const pointField& cellCentres,
171 const vectorField& faceCentres,
173
174 const scalarField& minOwnHeight,
175 const scalarField& minNeiHeight,
176
178) const
179{
180 // Clip correction vector if any pyramid becomes too small. Return number of
181 // cells clipped
182
183 const labelList& own = mesh_.faceOwner();
184 const labelList& nei = mesh_.faceNeighbour();
185
186 label nClipped = 0;
187 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
188 {
189 #ifdef WM_SPDP
190 const solveVector n(faceNormals[facei]);
191 const solveVector fc(faceCentres[facei]);
192 #else
193 const vector& n = faceNormals[facei];
194 const point& fc = faceCentres[facei];
195 #endif
196
197 const label ownCelli = own[facei];
198 if (correction[ownCelli] != vector::zero)
199 {
200 const solveVector ownCc(cellCentres[ownCelli]+correction[ownCelli]);
201 const scalar ownHeight = ((fc-ownCc)&n);
202 if (ownHeight < minOwnHeight[facei])
203 {
204 //Pout<< " internalface:" << fc
205 // << " own:" << ownCc
206 // << " pyrHeight:" << ownHeight
207 // << " minHeight:" << minOwnHeight[facei]
208 // << endl;
209 correction[ownCelli] = vector::zero;
210 nClipped++;
211 }
212 }
213
214 const label neiCelli = nei[facei];
215 if (correction[neiCelli] != vector::zero)
216 {
217 const solveVector neiCc(cellCentres[neiCelli]+correction[neiCelli]);
218 const scalar neiHeight = ((neiCc-fc)&n);
219 if (neiHeight < minNeiHeight[facei])
220 {
221 //Pout<< " internalface:" << fc
222 // << " nei:" << neiCc
223 // << " pyrHeight:" << neiHeight
224 // << " minHeight:" << minNeiHeight[facei]
225 // << endl;
226 correction[neiCelli] = vector::zero;
227 nClipped++;
228 }
229 }
230 }
231
232 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
233 {
234 #ifdef WM_SPDP
235 const solveVector n(faceNormals[facei]);
236 const solveVector fc(faceCentres[facei]);
237 #else
238 const vector& n = faceNormals[facei];
239 const point& fc = faceCentres[facei];
240 #endif
241
242 const label ownCelli = own[facei];
243 if (correction[ownCelli] != vector::zero)
244 {
245 const solveVector ownCc(cellCentres[ownCelli]+correction[ownCelli]);
246 const scalar ownHeight = ((fc-ownCc)&n);
247 if (ownHeight < minOwnHeight[facei])
248 {
249 //Pout<< " boundaryface:" << fc
250 // << " own:" << ownCc
251 // << " pyrHeight:" << ownHeight
252 // << " minHeight:" << minOwnHeight[facei]
253 // << endl;
254 correction[ownCelli] = vector::zero;
255 nClipped++;
256 }
257 }
258 }
259 return returnReduce(nClipped, sumOp<label>());
260}
261
262
265(
266 const pointField& cellCentres,
268 const scalarField& faceWeights
269) const
270{
271 const labelList& own = mesh_.faceOwner();
272 const labelList& nei = mesh_.faceNeighbour();
273
274
275 tmp<pointField> tcc(new pointField(mesh_.nCells(), Zero));
276 pointField& cc = tcc.ref();
277
278 Field<solveScalar> cellWeights(mesh_.nCells(), Zero);
279
280 // Internal faces
281 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
282 {
283 #ifdef WM_SPDP
284 const solveVector n(faceNormals[facei]);
285 #else
286 const vector& n = faceNormals[facei];
287 #endif
288 const point& ownCc = cellCentres[own[facei]];
289 const point& neiCc = cellCentres[nei[facei]];
290
291 solveVector d(neiCc-ownCc);
292
293 // 1. Normalise contribution. This increases actual non-ortho
294 // since it does not 'see' the tangential offset of neighbours
295 //neiCc = ownCc + (d&n)*n;
296
297 // 2. Remove normal contribution, i.e. get tangential vector
298 // (= non-ortho correction vector?)
299 d -= (d&n)*n;
300
301 // Apply half to both sides (as a correction)
302 // Note: should this be linear weights instead of 0.5?
303 const scalar w = 0.5*faceWeights[facei];
304 cc[own[facei]] += point(w*d);
305 cellWeights[own[facei]] += w;
306
307 cc[nei[facei]] -= point(w*d);
308 cellWeights[nei[facei]] += w;
309 }
310
311
312 // Boundary faces. Bypass stored cell centres
313 pointField neiCellCentres;
314 syncTools::swapBoundaryCellPositions(mesh_, cellCentres, neiCellCentres);
315
316 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
317 for (const polyPatch& pp : pbm)
318 {
319 if (pp.coupled())
320 {
321 const labelUList& fc = pp.faceCells();
322
323 forAll(fc, i)
324 {
325 const label meshFacei = pp.start()+i;
326 const label bFacei = meshFacei-mesh_.nInternalFaces();
327
328 #ifdef WM_SPDP
329 const solveVector n(faceNormals[meshFacei]);
330 #else
331 const vector& n = faceNormals[meshFacei];
332 #endif
333
334 const point& ownCc = cellCentres[fc[i]];
335 const point& neiCc = neiCellCentres[bFacei];
336
337 solveVector d(neiCc-ownCc);
338
339 // 1. Normalise contribution. This increases actual non-ortho
340 // since it does not 'see' the tangential offset of neighbours
341 //neiCc = ownCc + (d&n)*n;
342
343 // 2. Remove normal contribution, i.e. get tangential vector
344 // (= non-ortho correction vector?)
345 d -= (d&n)*n;
346
347 // Apply half to both sides (as a correction)
348 const scalar w = 0.5*faceWeights[meshFacei];
349 cc[fc[i]] += point(w*d);
350 cellWeights[fc[i]] += w;
351 }
352 }
353 }
354
355 // Now cc is still the correction vector. Add to cell original centres.
356 forAll(cc, celli)
357 {
358 if (cellWeights[celli] > VSMALL)
359 {
360 cc[celli] = cellCentres[celli] + cc[celli]/cellWeights[celli];
361 }
362 else
363 {
364 cc[celli] = cellCentres[celli];
365 }
366 }
367
368 return tcc;
369}
370
371
374(
375// const scalar ratio, // Amount of change in face-triangles area
376 const pointField& cellCentres,
377 const pointField& faceCentres,
379) const
380{
381 const labelList& own = mesh_.faceOwner();
382 const labelList& nei = mesh_.faceNeighbour();
383
384
385 tmp<pointField> tnewFc(new pointField(faceCentres));
386 pointField& newFc = tnewFc.ref();
387
388 // Internal faces
389 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
390 {
391 #ifdef WM_SPDP
392 const solveVector n(faceNormals[facei]);
393 const solveVector oldFc(faceCentres[facei]);
394 #else
395 const vector& n = faceNormals[facei];
396 const point& oldFc = faceCentres[facei];
397 #endif
398
399 const solveVector ownCc(cellCentres[own[facei]]);
400 const solveVector neiCc(cellCentres[nei[facei]]);
401
402 solveVector deltaCc(neiCc-ownCc);
403 solveVector deltaFc(oldFc-ownCc);
404
405 //solveVector d(neiCc-ownCc);
409 //
412 //d -= s*n;
413 //newFc[facei] = faceCentres[facei]+d;
414
415 // Get linear weight (normal distance to face)
416 const solveScalar f = (deltaFc&n)/(deltaCc&n);
417 const solveVector avgCc((1.0-f)*ownCc + f*neiCc);
418
419 solveVector d(avgCc-oldFc);
420 // Remove normal contribution, i.e. get tangential vector
421 // (= non-ortho correction vector?)
422 d -= (d&n)*n;
423
424// // Clip to limit change in
425// d *= ratio;
426
427
428 newFc[facei] = oldFc + d;
429 }
430
431
432 // Boundary faces. Bypass stored cell centres
433 pointField neiCellCentres;
434 syncTools::swapBoundaryCellPositions(mesh_, cellCentres, neiCellCentres);
435
436 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
437 for (const polyPatch& pp : pbm)
438 {
439 const labelUList& fc = pp.faceCells();
440
441 if (pp.coupled())
442 {
443 forAll(fc, i)
444 {
445 // Same as internal faces
446 const label facei = pp.start()+i;
447 const label bFacei = facei-mesh_.nInternalFaces();
448
449 #ifdef WM_SPDP
450 const solveVector n(faceNormals[facei]);
451 const solveVector oldFc(faceCentres[facei]);
452 #else
453 const vector& n = faceNormals[facei];
454 const point& oldFc = faceCentres[facei];
455 #endif
456
457 const solveVector ownCc(cellCentres[fc[i]]);
458 const solveVector neiCc(neiCellCentres[bFacei]);
459
460 solveVector deltaCc(neiCc-ownCc);
461 solveVector deltaFc(oldFc-ownCc);
462
463 // Get linear weight (normal distance to face)
464 const solveScalar f = (deltaFc&n)/(deltaCc&n);
465 const solveVector avgCc((1.0-f)*ownCc + f*neiCc);
466
467 solveVector d(avgCc-oldFc);
468 // Remove normal contribution, i.e. get tangential vector
469 // (= non-ortho correction vector?)
470 d -= (d&n)*n;
471
472 newFc[facei] = oldFc + d;
473 }
474 }
475 else
476 {
477 // Zero-grad?
478 forAll(fc, i)
479 {
480 const label facei = pp.start()+i;
481
482 #ifdef WM_SPDP
483 const solveVector n(faceNormals[facei]);
484 const solveVector oldFc(faceCentres[facei]);
485 #else
486 const vector& n = faceNormals[facei];
487 const point& oldFc = faceCentres[facei];
488 #endif
489
490 const solveVector ownCc(cellCentres[fc[i]]);
491
492 solveVector d(ownCc-oldFc);
493 // Remove normal contribution, i.e. get tangential vector
494 // (= non-ortho correction vector?)
495 d -= (d&n)*n;
496
497 newFc[facei] = oldFc+d;
498 }
499 }
500 }
501
502 return tnewFc;
503}
504
505
507(
508 const pointField& cellCentres,
510
511 scalarField& cosAngles,
512 scalarField& faceWeights
513) const
514{
515 cosAngles =
516 max
517 (
518 scalar(0),
519 min
520 (
521 scalar(1),
523 (
524 mesh_,
526 cellCentres
527 )
528 )
529 );
530
531
532 // Make weight: 0 for ortho faces, 1 for 90degrees non-ortho
533 //const scalarField faceWeights(scalar(1)-cosAngles);
534 faceWeights.setSize(cosAngles.size());
535 {
536 const scalar minCos = Foam::cos(degToRad(80));
537 const scalar maxCos = Foam::cos(degToRad(10));
538
539 forAll(cosAngles, facei)
540 {
541 const scalar cosAngle = cosAngles[facei];
542 if (cosAngle < minCos)
543 {
544 faceWeights[facei] = 1.0;
545 }
546 else if (cosAngle > maxCos)
547 {
548 faceWeights[facei] = 0.0;
549 }
550 else
551 {
552 faceWeights[facei] =
553 1.0-(cosAngle-minCos)/(maxCos-minCos);
554 }
555 }
556 }
557}
558
559
560// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
561
563(
564 const fvMesh& mesh,
565 const dictionary& dict
566)
567:
569 nIters_
570 (
571 dict.getCheckOrDefault<label>
572 (
573 "nIters",
574 1,
575 [](label val) { return val >= 0; }
576 )
577 ),
578 relax_
579 (
580 dict.getCheck<scalar>
581 (
582 "relax",
583 [](scalar val) { return val > 0 && val <= 1; }
584 )
585 ),
586 minRatio_
587 (
588 dict.getCheckOrDefault<scalar>
589 (
590 "minRatio",
591 0.5,
592 [](scalar val) { return val >= 0 && val <= 1; }
593 )
594 )
595{
596 if (debug)
597 {
598 Pout<< "averageNeighbourFvGeometryScheme :"
599 << " nIters:" << nIters_
600 << " relax:" << relax_
601 << " minRatio:" << minRatio_ << endl;
602 }
603
604 // Force local calculation
605 movePoints();
606}
607
608
609// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
610
612{
613 if (debug)
614 {
615 Pout<< "averageNeighbourFvGeometryScheme::movePoints() : "
616 << "recalculating primitiveMesh centres" << endl;
617 }
618
619 //if
620 //(
621 // !mesh_.hasCellCentres()
622 //&& !mesh_.hasFaceCentres()
623 //&& !mesh_.hasCellVolumes()
624 //&& !mesh_.hasFaceAreas()
625 //)
626 {
628
629 // Note: at this point the highAspectRatioFvGeometryScheme constructor
630 // will have already reset the primitive geometry!
631
632 vectorField faceAreas(mesh_.faceAreas());
633 const scalarField magFaceAreas(mag(faceAreas));
634 const vectorField faceNormals(faceAreas/magFaceAreas);
635
636
637 // Calculate aspectratio weights
638 // - 0 if aratio < minAspect_
639 // - 1 if aratio >= maxAspect_
640 scalarField cellWeight, faceWeight;
641 calcAspectRatioWeights(cellWeight, faceWeight);
642
643 // Relaxation
644 cellWeight *= relax_;
645 //faceWeight *= relax_;
646
647 // Calculate current pyramid heights
648 scalarField minOwnHeight;
649 scalarField minNeiHeight;
650 makePyrHeights
651 (
652 mesh_.cellCentres(),
653 mesh_.faceCentres(),
655
656 minOwnHeight,
657 minNeiHeight
658 );
659
660 // How much is the cell centre to vary inside the cell.
661 minOwnHeight *= minRatio_;
662 minNeiHeight *= minRatio_;
663
664
665
666 autoPtr<OBJstream> osPtr;
667 autoPtr<surfaceWriter> writerPtr;
668 if (debug)
669 {
670 osPtr.reset
671 (
672 new OBJstream
673 (
674 mesh_.time().timePath()
675 / "cellCentre_correction.obj"
676 )
677 );
678 Pout<< "averageNeighbourFvGeometryScheme::movePoints() :"
679 << " writing cell centre path to " << osPtr().name() << endl;
680
681
682 // Write current non-ortho
683 fileName outputDir
684 (
685 mesh_.time().globalPath()
687 / mesh_.pointsInstance()
688 );
689 outputDir.clean(); // Remove unneeded ".."
690 writerPtr = surfaceWriter::New
691 (
692 "ensight" //"vtk"
693 // options
694 );
695
696 // Use outputDir/TIME/surface-name
697 writerPtr->useTimeDir(true);
698
699 writerPtr->beginTime(mesh_.time());
700
701 writerPtr->open
702 (
703 mesh_.points(),
704 mesh_.faces(),
705 (outputDir / "cosAngle"),
706 true // merge parallel bits
707 );
708
709 writerPtr->endTime();
710 }
711
712
713 // Current cellCentres. These get adjusted to lower the
714 // non-orthogonality
715 pointField cellCentres(mesh_.cellCentres());
716
717 // Modify cell centres to be more in-line with the face normals
718 for (label iter = 0; iter < nIters_; iter++)
719 {
720 // Get neighbour average (weighted by face area). This gives
721 // preference to the dominant faces. However if the non-ortho
722 // is not caused by the dominant faces this moves to the wrong
723 // direction.
724 //tmp<pointField> tcc
725 //(
726 // averageNeighbourCentres
727 // (
728 // cellCentres,
729 // faceNormals,
730 // magFaceAreas
731 // )
732 //);
733
734 // Get neighbour average weighted by non-ortho. Question: how to
735 // weight boundary faces?
736 tmp<pointField> tcc;
737 {
738 scalarField cosAngles;
739 scalarField faceWeights;
740 makeNonOrthoWeights
741 (
742 cellCentres,
744
745 cosAngles,
746 faceWeights
747 );
748
749 if (writerPtr)
750 {
751 writerPtr->beginTime(instant(scalar(iter)));
752 writerPtr->write("cosAngles", cosAngles);
753 writerPtr->endTime();
754 }
755
756 if (debug)
757 {
758 forAll(cosAngles, facei)
759 {
760 if (cosAngles[facei] < Foam::cos(degToRad(85.0)))
761 {
762 Pout<< " face:" << facei
763 << " at:" << mesh_.faceCentres()[facei]
764 << " cosAngle:" << cosAngles[facei]
765 << " faceWeight:" << faceWeights[facei]
766 << endl;
767 }
768 }
769 }
770
771 tcc = averageNeighbourCentres
772 (
773 cellCentres,
775 faceWeights
776 );
777 }
778
779
780 // Calculate correction for cell centres. Leave low-aspect
781 // ratio cells unaffected (i.e. correction = 0)
782 vectorField correction(cellWeight*(tcc-cellCentres));
783
784 // Clip correction vector if pyramid becomes too small
785 const label nClipped = clipPyramids
786 (
787 cellCentres,
788 mesh_.faceCentres(),
790
791 minOwnHeight, // minimum owner pyramid height. Usually fraction
792 minNeiHeight, // of starting mesh
793
795 );
796
797 cellCentres += correction;
798
799 if (debug)
800 {
801 forAll(cellCentres, celli)
802 {
803 const point& cc = cellCentres[celli];
804 osPtr().write(linePointRef(cc-correction[celli], cc));
805 }
806
807 const scalarField magCorrection(mag(correction));
808 const scalarField faceOrthogonality
809 (
810 min
811 (
812 scalar(1),
814 (
815 mesh_,
816 faceAreas,
817 cellCentres
818 )
819 )
820 );
821 const scalarField nonOrthoAngle
822 (
824 (
825 Foam::acos(faceOrthogonality)
826 )
827 );
828 Pout<< " iter:" << iter
829 << " nClipped:" << nClipped
830 << " average displacement:" << gAverage(magCorrection)
831 << " non-ortho angle : average:" << gAverage(nonOrthoAngle)
832 << " max:" << gMax(nonOrthoAngle) << endl;
833 }
834 }
835
837 (
838 averageCentres
839 (
840 cellCentres,
841 mesh_.faceCentres(),
843 )
844 );
845 vectorField faceCorrection(faceWeight*(tfc-mesh_.faceCentres()));
846 // Clip correction vector to not have min triangle shrink
847 // by more than minRatio
848 clipFaceTet
849 (
850 minRatio_,
851 mesh_.faceCentres(),
853 faceCorrection
854 );
855 vectorField faceCentres(mesh_.faceCentres() + faceCorrection);
856
857 if (debug)
858 {
859 Pout<< "averageNeighbourFvGeometryScheme::movePoints() :"
860 << " averageNeighbour weight"
861 << " max:" << gMax(cellWeight) << " min:" << gMin(cellWeight)
862 << " average:" << gAverage(cellWeight) << endl;
863
864 // Dump lines from old to new location
865 const fileName tp(mesh_.time().timePath());
866 mkDir(tp);
867 OBJstream str(tp/"averageNeighbourCellCentres.obj");
868 Pout<< "Writing lines from old to new cell centre to " << str.name()
869 << endl;
870 forAll(mesh_.cellCentres(), celli)
871 {
872 const point& oldCc = mesh_.cellCentres()[celli];
873 const point& newCc = cellCentres[celli];
874 str.write(linePointRef(oldCc, newCc));
875 }
876 }
877 if (debug)
878 {
879 // Dump lines from old to new location
880 const fileName tp(mesh_.time().timePath());
881 OBJstream str(tp/"averageFaceCentres.obj");
882 Pout<< "Writing lines from old to new face centre to " << str.name()
883 << endl;
884 forAll(mesh_.faceCentres(), facei)
885 {
886 const point& oldFc = mesh_.faceCentres()[facei];
887 const point& newFc = faceCentres[facei];
888 str.write(linePointRef(oldFc, newFc));
889 }
890 }
891
892 scalarField cellVolumes(mesh_.cellVolumes());
893
894 // Store on primitiveMesh
895 //const_cast<fvMesh&>(mesh_).clearGeom();
896 const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
897 (
898 std::move(faceCentres),
899 std::move(faceAreas),
900 std::move(cellCentres),
901 std::move(cellVolumes)
902 );
903 }
904}
905
906
907// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
OFstream that keeps track of vertices.
Definition: OBJstream.H:61
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
const T & fcValue(const label i) const
Return forward circular value (ie, next value in the list)
Definition: UListI.H:74
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Geometry calculation scheme to minimise non-orthogonality/.
label clipFaceTet(const scalar minRatio, const vectorField &faceCentres, const vectorField &faceNormals, vectorField &faceCorrection) const
tmp< pointField > averageCentres(const pointField &cellCentres, const pointField &faceCentres, const vectorField &faceNormals) const
virtual void movePoints()
Do what is necessary if the mesh has moved.
void makePyrHeights(const pointField &cellCentres, const vectorField &faceCentres, const vectorField &faceNormals, scalarField &ownHeight, scalarField &neiHeight) const
Calculate pyramid heights.
tmp< pointField > averageNeighbourCentres(const pointField &cellCentres, const vectorField &faceNormals, const scalarField &faceWeights) const
Average neighbouring cell centres to minimise non-ortho.
label clipPyramids(const pointField &cellCentres, const vectorField &faceCentres, const vectorField &faceNormals, const scalarField &minOwnHeight, const scalarField &minNeiHeight, vectorField &correction) const
void makeNonOrthoWeights(const pointField &cellCentres, const vectorField &faceNormals, scalarField &cosAngles, scalarField &faceWeights) const
Make weights based on non-orthogonality.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getCheck(const word &keyword, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A class for handling file names.
Definition: fileName.H:76
static bool clean(std::string &str)
Definition: fileName.C:199
static word outputPrefix
Directory prefix.
Abstract base class for geometry calculation schemes.
const fvMesh & mesh_
Hold reference to mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Geometry calculation scheme with automatic stabilisation for high-aspect ratio cells.
virtual void movePoints()
Do what is necessary if the mesh has moved.
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition: instant.H:56
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Definition: polyMeshTools.C:37
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label nFaces() const noexcept
Number of mesh faces.
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
static void swapBoundaryCellPositions(const polyMesh &mesh, const UList< point > &cellData, List< point > &neighbourCellData)
Swap to obtain neighbour cell positions for all boundary faces.
Definition: syncTools.C:34
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
surfaceVectorField faceNormals(mesh.Sf()/mesh.magSf())
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:515
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Type gMin(const FieldField< Field, Type > &f)
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Type gMax(const FieldField< Field, Type > &f)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59
Unit conversion functions.