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