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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
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 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(averageNeighbourFvGeometryScheme, 0);
44  (
45  fvGeometryScheme,
46  averageNeighbourFvGeometryScheme,
47  dict
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
55 (
56  const scalar minRatio,
57  const vectorField& faceCentres,
58  const vectorField& faceNormals,
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,
134  const vectorField& faceNormals,
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]];
163  ownHeight[facei] = ((fc-ownCc)&n);
164  }
165 }
166 
167 
169 (
170  const pointField& cellCentres,
171  const vectorField& faceCentres,
172  const vectorField& faceNormals,
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,
267  const vectorField& faceNormals,
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  }
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,
378  const vectorField& faceNormals
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,
509  const vectorField& faceNormals,
510 
511  scalarField& cosAngles,
512  scalarField& faceWeights
513 ) const
514 {
515  cosAngles = clamp
516  (
517  polyMeshTools::faceOrthogonality(mesh_, faceNormals, cellCentres),
518  zero_one{}
519  );
520 
521  // Make weight: 0 for ortho faces, 1 for 90degrees non-ortho
522  //const scalarField faceWeights(scalar(1)-cosAngles);
523  faceWeights.setSize(cosAngles.size());
524  {
525  const scalar minCos = Foam::cos(degToRad(80));
526  const scalar maxCos = Foam::cos(degToRad(10));
527 
528  forAll(cosAngles, facei)
529  {
530  const scalar cosAngle = cosAngles[facei];
531  if (cosAngle < minCos)
532  {
533  faceWeights[facei] = 1.0;
534  }
535  else if (cosAngle > maxCos)
536  {
537  faceWeights[facei] = 0.0;
538  }
539  else
540  {
541  faceWeights[facei] =
542  1.0-(cosAngle-minCos)/(maxCos-minCos);
543  }
544  }
545  }
546 }
547 
548 
549 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
550 
551 Foam::averageNeighbourFvGeometryScheme::averageNeighbourFvGeometryScheme
552 (
553  const fvMesh& mesh,
554  const dictionary& dict
555 )
556 :
557  highAspectRatioFvGeometryScheme(mesh, dict),
558  nIters_
559  (
560  dict.getCheckOrDefault<label>
561  (
562  "nIters",
563  1,
564  [](label val) { return val >= 0; }
565  )
566  ),
567  relax_
568  (
569  dict.getCheck<scalar>
570  (
571  "relax",
572  [](scalar val) { return val > 0 && val <= 1; }
573  )
574  ),
575  minRatio_
576  (
577  dict.getCheckOrDefault<scalar>
578  (
579  "minRatio",
580  0.5,
581  [](scalar val) { return val >= 0 && val <= 1; }
582  )
583  )
584 {
585  if (debug)
586  {
587  Pout<< "averageNeighbourFvGeometryScheme :"
588  << " nIters:" << nIters_
589  << " relax:" << relax_
590  << " minRatio:" << minRatio_ << endl;
591  }
592 
593  // Force local calculation
594  movePoints();
595 }
596 
597 
598 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
599 
601 {
602  if (debug)
603  {
604  Pout<< "averageNeighbourFvGeometryScheme::movePoints() : "
605  << "recalculating primitiveMesh centres" << endl;
606  }
607 
608  //if
609  //(
610  // !mesh_.hasCellCentres()
611  //&& !mesh_.hasFaceCentres()
612  //&& !mesh_.hasCellVolumes()
613  //&& !mesh_.hasFaceAreas()
614  //)
615  {
617 
618  // Note: at this point the highAspectRatioFvGeometryScheme constructor
619  // will have already reset the primitive geometry!
620 
621  vectorField faceAreas(mesh_.faceAreas());
622  const scalarField magFaceAreas(mag(faceAreas));
623  const vectorField faceNormals(faceAreas/magFaceAreas);
624 
625 
626  // Calculate aspectratio weights
627  // - 0 if aratio < minAspect_
628  // - 1 if aratio >= maxAspect_
629  scalarField cellWeight, faceWeight;
630  calcAspectRatioWeights(cellWeight, faceWeight);
631 
632  // Relaxation
633  cellWeight *= relax_;
634  //faceWeight *= relax_;
635 
636  // Calculate current pyramid heights
637  scalarField minOwnHeight;
638  scalarField minNeiHeight;
639  makePyrHeights
640  (
641  mesh_.cellCentres(),
642  mesh_.faceCentres(),
643  faceNormals,
644 
645  minOwnHeight,
646  minNeiHeight
647  );
648 
649  // How much is the cell centre to vary inside the cell.
650  minOwnHeight *= minRatio_;
651  minNeiHeight *= minRatio_;
652 
653 
654 
655  autoPtr<OBJstream> osPtr;
656  autoPtr<surfaceWriter> writerPtr;
657  if (debug)
658  {
659  osPtr.reset
660  (
661  new OBJstream
662  (
663  mesh_.time().timePath()
664  / "cellCentre_correction.obj"
665  )
666  );
667  Pout<< "averageNeighbourFvGeometryScheme::movePoints() :"
668  << " writing cell centre path to " << osPtr().name() << endl;
669 
670 
671  // Write current non-ortho
672  fileName outputDir
673  (
674  mesh_.time().globalPath()
676  / mesh_.pointsInstance()
677  );
678  outputDir.clean(); // Remove unneeded ".."
679  writerPtr = surfaceWriter::New
680  (
681  "ensight" //"vtk"
682  // options
683  );
684 
685  // Use outputDir/TIME/surface-name
686  writerPtr->useTimeDir(true);
687 
688  writerPtr->beginTime(mesh_.time());
689 
690  writerPtr->open
691  (
692  mesh_.points(),
693  mesh_.faces(),
694  (outputDir / "cosAngle"),
695  true // merge parallel bits
696  );
697 
698  writerPtr->endTime();
699  }
700 
701 
702  // Current cellCentres. These get adjusted to lower the
703  // non-orthogonality
704  pointField cellCentres(mesh_.cellCentres());
705 
706  // Modify cell centres to be more in-line with the face normals
707  for (label iter = 0; iter < nIters_; iter++)
708  {
709  // Get neighbour average (weighted by face area). This gives
710  // preference to the dominant faces. However if the non-ortho
711  // is not caused by the dominant faces this moves to the wrong
712  // direction.
713  //tmp<pointField> tcc
714  //(
715  // averageNeighbourCentres
716  // (
717  // cellCentres,
718  // faceNormals,
719  // magFaceAreas
720  // )
721  //);
722 
723  // Get neighbour average weighted by non-ortho. Question: how to
724  // weight boundary faces?
725  tmp<pointField> tcc;
726  {
727  scalarField cosAngles;
728  scalarField faceWeights;
729  makeNonOrthoWeights
730  (
731  cellCentres,
732  faceNormals,
733 
734  cosAngles,
735  faceWeights
736  );
737 
738  if (writerPtr)
739  {
740  writerPtr->beginTime(instant(scalar(iter)));
741  writerPtr->write("cosAngles", cosAngles);
742  writerPtr->endTime();
743  }
744 
745  if (debug)
746  {
747  forAll(cosAngles, facei)
748  {
749  if (cosAngles[facei] < Foam::cos(degToRad(85.0)))
750  {
751  Pout<< " face:" << facei
752  << " at:" << mesh_.faceCentres()[facei]
753  << " cosAngle:" << cosAngles[facei]
754  << " faceWeight:" << faceWeights[facei]
755  << endl;
756  }
757  }
758  }
759 
760  tcc = averageNeighbourCentres
761  (
762  cellCentres,
763  faceNormals,
764  faceWeights
765  );
766  }
767 
768 
769  // Calculate correction for cell centres. Leave low-aspect
770  // ratio cells unaffected (i.e. correction = 0)
771  vectorField correction(cellWeight*(tcc-cellCentres));
772 
773  // Clip correction vector if pyramid becomes too small
774  const label nClipped = clipPyramids
775  (
776  cellCentres,
777  mesh_.faceCentres(),
778  faceNormals,
779 
780  minOwnHeight, // minimum owner pyramid height. Usually fraction
781  minNeiHeight, // of starting mesh
782 
783  correction
784  );
785 
786  cellCentres += correction;
787 
788  if (debug)
789  {
790  forAll(cellCentres, celli)
791  {
792  const point& cc = cellCentres[celli];
793  osPtr().writeLine(cc-correction[celli], cc);
794  }
795 
796  const scalarField magCorrection(mag(correction));
797  const scalarField faceOrthogonality
798  (
799  min
800  (
801  scalar(1),
803  (
804  mesh_,
805  faceAreas,
806  cellCentres
807  )
808  )
809  );
810  const scalarField nonOrthoAngle
811  (
812  radToDeg
813  (
814  Foam::acos(faceOrthogonality)
815  )
816  );
817  Pout<< " iter:" << iter
818  << " nClipped:" << nClipped
819  << " average displacement:" << gAverage(magCorrection)
820  << " non-ortho angle : average:" << gAverage(nonOrthoAngle)
821  << " max:" << gMax(nonOrthoAngle) << endl;
822  }
823  }
824 
825  tmp<pointField> tfc
826  (
827  averageCentres
828  (
829  cellCentres,
830  mesh_.faceCentres(),
831  faceNormals
832  )
833  );
834  vectorField faceCorrection(faceWeight*(tfc-mesh_.faceCentres()));
835  // Clip correction vector to not have min triangle shrink
836  // by more than minRatio
837  clipFaceTet
838  (
839  minRatio_,
840  mesh_.faceCentres(),
841  faceNormals,
842  faceCorrection
843  );
844  vectorField faceCentres(mesh_.faceCentres() + faceCorrection);
845 
846  if (debug)
847  {
848  Pout<< "averageNeighbourFvGeometryScheme::movePoints() :"
849  << " averageNeighbour weight"
850  << " max:" << gMax(cellWeight) << " min:" << gMin(cellWeight)
851  << " average:" << gAverage(cellWeight) << endl;
852 
853  // Dump lines from old to new location
854  const fileName tp(mesh_.time().timePath());
855  mkDir(tp);
856  OBJstream str(tp/"averageNeighbourCellCentres.obj");
857  Pout<< "Writing lines from old to new cell centre to " << str.name()
858  << endl;
859  forAll(mesh_.cellCentres(), celli)
860  {
861  const point& oldCc = mesh_.cellCentres()[celli];
862  const point& newCc = cellCentres[celli];
863  str.writeLine(oldCc, newCc);
864  }
865  }
866  if (debug)
867  {
868  // Dump lines from old to new location
869  const fileName tp(mesh_.time().timePath());
870  OBJstream str(tp/"averageFaceCentres.obj");
871  Pout<< "Writing lines from old to new face centre to " << str.name()
872  << endl;
873  forAll(mesh_.faceCentres(), facei)
874  {
875  const point& oldFc = mesh_.faceCentres()[facei];
876  const point& newFc = faceCentres[facei];
877  str.writeLine(oldFc, newFc);
878  }
879  }
880 
881  scalarField cellVolumes(mesh_.cellVolumes());
882 
883  // Store on primitiveMesh
884  //const_cast<fvMesh&>(mesh_).clearGeom();
885  const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
886  (
887  std::move(faceCentres),
888  std::move(faceAreas),
889  std::move(cellCentres),
890  std::move(cellVolumes)
891  );
892  }
893 }
894 
895 
896 // ************************************************************************* //
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
const polyBoundaryMesh & pbm
dictionary dict
dimensionedScalar acos(const dimensionedScalar &ds)
const fvMesh & mesh_
Hold reference to mesh.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type gMin(const FieldField< Field, Type > &f)
tmp< pointField > averageCentres(const pointField &cellCentres, const pointField &faceCentres, const vectorField &faceNormals) const
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
Unit conversion functions.
T getCheck(const word &keyword, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T with additional checking FatalIOError if not found, or if the number of tokens is...
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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:27
label nFaces() const noexcept
Number of mesh faces.
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Vector< solveScalar > solveVector
Definition: vector.H:60
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
virtual void movePoints()
Do what is necessary if the mesh has moved.
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
dimensionedScalar cos(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
tmp< pointField > averageNeighbourCentres(const pointField &cellCentres, const vectorField &faceNormals, const scalarField &faceWeights) const
Average neighbouring cell centres to minimise non-ortho.
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross-product operators.
Definition: Vector.H:58
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
label clipFaceTet(const scalar minRatio, const vectorField &faceCentres, const vectorField &faceNormals, vectorField &faceCorrection) const
Clip face-centre correction vector if new triangle area would get below min. Return number of clipped...
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
void makeNonOrthoWeights(const pointField &cellCentres, const vectorField &faceNormals, scalarField &cosAngles, scalarField &faceWeights) const
Make weights based on non-orthogonality.
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate orthogonality field. (1 for fully orthogonal, < 1 for non-orthogonal)
Definition: polyMeshTools.C:30
vector point
Point is a vector.
Definition: point.H:37
Type gAverage(const FieldField< Field, Type > &f)
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
label clipPyramids(const pointField &cellCentres, const vectorField &faceCentres, const vectorField &faceNormals, const scalarField &minOwnHeight, const scalarField &minNeiHeight, vectorField &correction) const
Clip correction vector if new pyramid height would get below min. Return number of clipped cells...
static word outputPrefix
Directory prefix.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
virtual void movePoints()
Do what is necessary if the mesh has moved.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:80
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
void makePyrHeights(const pointField &cellCentres, const vectorField &faceCentres, const vectorField &faceNormals, scalarField &ownHeight, scalarField &neiHeight) const
Calculate pyramid heights.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127