highAspectRatioFvGeometryScheme.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 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 "triangle.H"
32 #include "syncTools.H"
33 #include "cellAspectRatio.H"
34 #include "emptyPolyPatch.H"
35 #include "wedgePolyPatch.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(highAspectRatioFvGeometryScheme, 0);
43  (
44  fvGeometryScheme,
45  highAspectRatioFvGeometryScheme,
46  dict
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 //void Foam::highAspectRatioFvGeometryScheme::cellClosedness
54 //(
55 // const vectorField& areas,
56 // const scalarField& vols,
57 // const tensorField& cellCoords,
58 //
59 // scalarField& aratio
60 //) const
61 //{
62 // // From primitiveMeshTools::cellClosedness:
63 // // calculate aspect ratio in given direction
64 // const labelList& own = mesh_.faceOwner();
65 // const labelList& nei = mesh_.faceNeighbour();
66 //
67 // // Loop through cell faces and sum up the face area vectors for each cell.
68 // // This should be zero in all vector components
69 //
70 // vectorField sumMagClosed(mesh_.nCells(), Zero);
71 //
72 // forAll(own, facei)
73 // {
74 // // Add to owner
75 // vector& v = sumMagClosed[own[facei]];
76 // v += cmptMag(cellCoords[own[facei]] & areas[facei]);
77 // }
78 //
79 // forAll(nei, facei)
80 // {
81 // // Subtract from neighbour
82 // vector& v = sumMagClosed[nei[facei]];
83 // v += cmptMag(cellCoords[nei[facei]] & areas[facei]);
84 // }
85 //
86 // // Check the sums
87 // aratio.setSize(mesh_.nCells());
88 //
89 // forAll(sumMagClosed, celli)
90 // {
91 // // Calculate the aspect ration as the maximum of Cartesian component
92 // // aspect ratio to the total area hydraulic area aspect ratio
93 // scalar minCmpt = VGREAT;
94 // scalar maxCmpt = -VGREAT;
95 // for (direction dir = 0; dir < vector::nComponents; dir++)
96 // {
97 // minCmpt = min(minCmpt, sumMagClosed[celli][dir]);
98 // maxCmpt = max(maxCmpt, sumMagClosed[celli][dir]);
99 // }
100 //
101 // scalar aspectRatio = maxCmpt/(minCmpt + ROOTVSMALL);
102 // const scalar v = max(ROOTVSMALL, vols[celli]);
103 //
104 // aspectRatio = max
105 // (
106 // aspectRatio,
107 // 1.0/6.0*cmptSum(sumMagClosed[celli])/Foam::pow(v, 2.0/3.0)
108 // );
109 //
110 // aratio[celli] = aspectRatio;
111 // }
112 //}
113 //
114 //
115 //void Foam::highAspectRatioFvGeometryScheme::cellDirections
116 //(
117 // tensorField& T,
118 // vectorField& lambda
119 //) const
120 //{
121 // // Calculate principal directions in increasing order
122 //
123 // T.setSize(mesh_.nCells());
124 // lambda.setSize(mesh_.nCells());
125 //
126 // forAll(T, celli)
127 // {
128 // tensor J = Zero;
129 // {
130 // const List<tetIndices> cellTets
131 // (
132 // polyMeshTetDecomposition::cellTetIndices
133 // (
134 // mesh_,
135 // celli
136 // )
137 // );
138 // triFaceList faces(cellTets.size());
139 // forAll(cellTets, cTI)
140 // {
141 // faces[cTI] = cellTets[cTI].faceTriIs(mesh_);
142 // }
143 //
144 // scalar m = 0.0;
145 // vector cM = Zero;
146 // J = Zero;
147 // momentOfInertia::massPropertiesShell
148 // (
149 // mesh_.points(),
150 // faces,
151 // 1.0,
152 // m,
153 // cM,
154 // J
155 // );
156 // }
157 //
158 // lambda[celli] = Foam::eigenValues(J);
159 // T[celli] = Foam::eigenVectors(J, lambda[celli]);
160 // }
161 //}
162 
164 (
165  scalarField& cellWeight,
166  scalarField& faceWeight
167 ) const
168 {
169  //scalarField aratio;
170  //{
171  // tensorField principalDirections;
172  // vectorField lambdas;
173  // cellDirections(principalDirections, lambdas);
174  //
175  // cellClosedness
176  // (
177  // mesh_.faceAreas(),
178  // mesh_.cellVolumes(),
179  // principalDirections,
180  // aratio
181  // );
182  //}
183  const cellAspectRatio aratio(mesh_);
184 
185  // Weighting for correction
186  // - 0 if aratio < minAspect_
187  // - 1 if aratio >= maxAspect_
188 
189  scalar delta(maxAspect_-minAspect_);
190  if (delta < ROOTVSMALL)
191  {
192  delta = SMALL;
193  }
194 
195  cellWeight = clamp((aratio-minAspect_)/delta, zero_one{});
196 
197  faceWeight.setSize(mesh_.nFaces());
198 
199  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
200  {
201  const label own = mesh_.faceOwner()[facei];
202  const label nei = mesh_.faceNeighbour()[facei];
203  faceWeight[facei] = max(cellWeight[own], cellWeight[nei]);
204  }
205  scalarField nbrCellWeight;
207  (
208  mesh_,
209  cellWeight,
210  nbrCellWeight
211  );
212  for
213  (
214  label facei = mesh_.nInternalFaces();
215  facei < mesh_.nFaces();
216  facei++
217  )
218  {
219  const label own = mesh_.faceOwner()[facei];
220  const label bFacei = facei-mesh_.nInternalFaces();
221  faceWeight[facei] = max(cellWeight[own], nbrCellWeight[bFacei]);
222  }
223 }
224 
225 
227 (
228  const polyMesh& mesh,
229  const pointField& p,
230  const pointField& faceAreas,
231  const scalarField& magFaceAreas,
232  pointField& faceCentres,
233  pointField& cellCentres
234 )
235 {
236  if (debug)
237  {
238  Pout<< "highAspectRatioFvGeometryScheme::makeAverageCentres() : "
239  << "calculating weighted average face/cell centre" << endl;
240  }
241 
242  const faceList& fs = mesh.faces();
243 
244  // Start off from primitiveMesh faceCentres (preserved for triangles)
245  faceCentres.setSize(mesh.nFaces());
246 
247  forAll(fs, facei)
248  {
249  const labelList& f = fs[facei];
250  const label nPoints = f.size();
251 
252  if (nPoints == 3)
253  {
254  faceCentres[facei] = triPointRef::centre(p[f[0]], p[f[1]], p[f[2]]);
255  }
256  else
257  {
258  solveScalar sumA = 0.0;
259  solveVector sumAc = Zero;
260 
261  for (label pi = 0; pi < nPoints; pi++)
262  {
263  const label nextPi(pi == nPoints-1 ? 0 : pi+1);
264  const solveVector nextPoint(p[f[nextPi]]);
265  const solveVector thisPoint(p[f[pi]]);
266 
267  const solveVector eMid = 0.5*(thisPoint+nextPoint);
268  const solveScalar a = mag(nextPoint-thisPoint);
269 
270  sumAc += a*eMid;
271  sumA += a;
272  }
273  // This is to deal with zero-area faces. Mark very small faces
274  // to be detected in e.g. processorPolyPatch.
275  if (sumA >= ROOTVSMALL)
276  {
277  faceCentres[facei] = sumAc/sumA;
278  }
279  else
280  {
281  // Unweighted average of points
282  sumAc = Zero;
283  for (label pi = 0; pi < nPoints; pi++)
284  {
285  sumAc += static_cast<solveVector>(p[f[pi]]);
286  }
287  faceCentres[facei] = sumAc/nPoints;
288  }
289  }
290  }
291 
292 
293  cellCentres.setSize(mesh.nCells());
294  cellCentres = Zero;
295  {
296  const labelList& own = mesh.faceOwner();
297  const labelList& nei = mesh.faceNeighbour();
298 
299  Field<solveScalar> cellWeights(mesh.nCells(), Zero);
300  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
301  {
302  const solveScalar magfA(magFaceAreas[facei]);
303  const vector weightedFc(magfA*faceCentres[facei]);
304 
305  // Accumulate area-weighted face-centre
306  cellCentres[own[facei]] += weightedFc;
307  cellCentres[nei[facei]] += weightedFc;
308 
309  // Accumulate weights
310  cellWeights[own[facei]] += magfA;
311  cellWeights[nei[facei]] += magfA;
312  }
313 
314  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
315  for (const polyPatch& pp : pbm)
316  {
317  if (!isA<emptyPolyPatch>(pp) && !isA<wedgePolyPatch>(pp))
318  {
319  for
320  (
321  label facei = pp.start();
322  facei < pp.start()+pp.size();
323  facei++
324  )
325  {
326  const solveScalar magfA(magFaceAreas[facei]);
327  const vector weightedFc(magfA*faceCentres[facei]);
328 
329  cellCentres[own[facei]] += weightedFc;
330  cellWeights[own[facei]] += magfA;
331  }
332  }
333  }
334 
335  forAll(cellCentres, celli)
336  {
337  if (mag(cellWeights[celli]) > VSMALL)
338  {
339  cellCentres[celli] /= cellWeights[celli];
340  }
341  }
342  }
343 }
344 
345 
346 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
347 
348 Foam::highAspectRatioFvGeometryScheme::highAspectRatioFvGeometryScheme
349 (
350  const fvMesh& mesh,
351  const dictionary& dict
352 )
353 :
354  basicFvGeometryScheme(mesh, dict),
355  minAspect_(dict.get<scalar>("minAspect")),
356  maxAspect_(dict.get<scalar>("maxAspect"))
357 {
358  if (maxAspect_ < minAspect_)
359  {
361  << "minAspect " << minAspect_
362  << " has to be less than maxAspect " << maxAspect_
363  << exit(FatalIOError);
364  }
365  if (minAspect_ < 0 || maxAspect_ < 0)
366  {
368  << "Illegal aspect ratio : minAspect:" << minAspect_
369  << " maxAspect:" << maxAspect_
370  << exit(FatalIOError);
371  }
372 
373  // Force local calculation
374  movePoints();
375 }
376 
377 
378 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
379 
381 {
382  //basicFvGeometryScheme::movePoints();
384 
385  if (debug)
386  {
387  Pout<< "highAspectRatioFvGeometryScheme::movePoints() : "
388  << "recalculating primitiveMesh centres" << endl;
389  }
390 
391  if
392  (
393  !mesh_.hasCellCentres()
394  && !mesh_.hasFaceCentres()
395  && !mesh_.hasCellVolumes()
396  && !mesh_.hasFaceAreas()
397  )
398  {
399  // Use lower level to calculate the geometry
400  const_cast<fvMesh&>(mesh_).primitiveMesh::updateGeom();
401 
402  pointField avgFaceCentres;
403  pointField avgCellCentres;
404  makeAverageCentres
405  (
406  mesh_,
407  mesh_.points(),
408  mesh_.faceAreas(),
409  mag(mesh_.faceAreas()),
410  avgFaceCentres,
411  avgCellCentres
412  );
413 
414 
415  // Calculate aspectratio weights
416  // - 0 if aratio < minAspect_
417  // - 1 if aratio >= maxAspect_
418  scalarField cellWeight, faceWeight;
419  calcAspectRatioWeights(cellWeight, faceWeight);
420 
421 
422  // Weight with average ones
423  vectorField faceCentres
424  (
425  (1.0-faceWeight)*mesh_.faceCentres()
426  + faceWeight*avgFaceCentres
427  );
428  vectorField cellCentres
429  (
430  (1.0-cellWeight)*mesh_.cellCentres()
431  + cellWeight*avgCellCentres
432  );
433 
434 
435  if (debug)
436  {
437  Pout<< "highAspectRatioFvGeometryScheme::movePoints() :"
438  << " highAspectRatio weight"
439  << " max:" << gMax(cellWeight) << " min:" << gMin(cellWeight)
440  << " average:" << gAverage(cellWeight) << endl;
441  }
442 
443  vectorField faceAreas(mesh_.faceAreas());
444  scalarField cellVolumes(mesh_.cellVolumes());
445 
446  // Store on primitiveMesh
447  //const_cast<fvMesh&>(mesh_).clearGeom();
448  const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
449  (
450  std::move(faceCentres),
451  std::move(faceAreas),
452  std::move(cellCentres),
453  std::move(cellVolumes)
454  );
455  }
456 }
457 
458 
459 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
scalar delta
const polyBoundaryMesh & pbm
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const fvMesh & mesh_
Hold reference to mesh.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type gMin(const FieldField< Field, Type > &f)
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
(Rough approximation of) cell aspect ratio
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void makeAverageCentres(const polyMesh &mesh, const pointField &points, const pointField &faceAreas, const scalarField &magFaceAreas, pointField &faceCentres, pointField &cellCentres)
Helper : calculate (weighted) average face and cell centres.
label nFaces() const noexcept
Number of mesh faces.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Vector< solveScalar > solveVector
Definition: vector.H:60
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
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.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
Point centre() const
Return centre (centroid)
Definition: triangleI.H:126
virtual void updateGeom()
Update all geometric data.
constexpr scalar pi(M_PI)
label nInternalFaces() const noexcept
Number of internal faces.
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:65
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
void calcAspectRatioWeights(scalarField &cellWeight, scalarField &faceWeight) const
Calculate cell and face weight. Is 0 for cell < minAspect, 1 for.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
label nCells() const noexcept
Number of mesh cells.
Type gAverage(const FieldField< Field, Type > &f)
virtual void movePoints()
Update basic geometric properties from provided points.
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127