polyMeshTools.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
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 
28 #include "polyMeshTools.H"
29 #include "syncTools.H"
30 #include "pyramidPointFaceRef.H"
31 #include "primitiveMeshTools.H"
32 
33 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34 
36 (
37  const polyMesh& mesh,
38  const vectorField& areas,
39  const vectorField& cc
40 )
41 {
42  const labelList& own = mesh.faceOwner();
43  const labelList& nei = mesh.faceNeighbour();
44  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
45 
46  tmp<scalarField> tortho(new scalarField(mesh.nFaces(), 1.0));
47  scalarField& ortho = tortho.ref();
48 
49  // Internal faces
50  forAll(nei, facei)
51  {
52  ortho[facei] = primitiveMeshTools::faceOrthogonality
53  (
54  cc[own[facei]],
55  cc[nei[facei]],
56  areas[facei]
57  );
58  }
59 
60 
61  // Coupled faces
62 
63  pointField neighbourCc;
64  syncTools::swapBoundaryCellPositions(mesh, cc, neighbourCc);
65 
66  forAll(pbm, patchi)
67  {
68  const polyPatch& pp = pbm[patchi];
69  if (pp.coupled())
70  {
71  forAll(pp, i)
72  {
73  label facei = pp.start() + i;
74  label bFacei = facei - mesh.nInternalFaces();
75 
76  ortho[facei] = primitiveMeshTools::faceOrthogonality
77  (
78  cc[own[facei]],
79  neighbourCc[bFacei],
80  areas[facei]
81  );
82  }
83  }
84  }
85 
86  return tortho;
87 }
88 
89 
91 (
92  const polyMesh& mesh,
93  const pointField& p,
94  const vectorField& fCtrs,
95  const vectorField& fAreas,
96  const vectorField& cellCtrs
97 )
98 {
99  const labelList& own = mesh.faceOwner();
100  const labelList& nei = mesh.faceNeighbour();
101  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
102 
103  tmp<scalarField> tskew(new scalarField(mesh.nFaces()));
104  scalarField& skew = tskew.ref();
105 
106  forAll(nei, facei)
107  {
108  skew[facei] = primitiveMeshTools::faceSkewness
109  (
110  mesh,
111  p,
112  fCtrs,
113  fAreas,
114 
115  facei,
116  cellCtrs[own[facei]],
117  cellCtrs[nei[facei]]
118  );
119  }
120 
121 
122  // Boundary faces: consider them to have only skewness error.
123  // (i.e. treat as if mirror cell on other side)
124 
125  pointField neighbourCc;
126  syncTools::swapBoundaryCellPositions(mesh, cellCtrs, neighbourCc);
127 
128  forAll(pbm, patchi)
129  {
130  const polyPatch& pp = pbm[patchi];
131  if (pp.coupled())
132  {
133  forAll(pp, i)
134  {
135  label facei = pp.start() + i;
136  label bFacei = facei - mesh.nInternalFaces();
137 
138  skew[facei] = primitiveMeshTools::faceSkewness
139  (
140  mesh,
141  p,
142  fCtrs,
143  fAreas,
144 
145  facei,
146  cellCtrs[own[facei]],
147  neighbourCc[bFacei]
148  );
149  }
150  }
151  else
152  {
153  forAll(pp, i)
154  {
155  label facei = pp.start() + i;
156 
157  skew[facei] = primitiveMeshTools::boundaryFaceSkewness
158  (
159  mesh,
160  p,
161  fCtrs,
162  fAreas,
163 
164  facei,
165  cellCtrs[own[facei]]
166  );
167  }
168  }
169  }
170 
171  return tskew;
172 }
173 
174 
176 (
177  const polyMesh& mesh,
178  const vectorField& fCtrs,
179  const vectorField& fAreas,
180  const vectorField& cellCtrs
181 )
182 {
183  const labelList& own = mesh.faceOwner();
184  const labelList& nei = mesh.faceNeighbour();
185  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
186 
187  tmp<scalarField> tweight(new scalarField(mesh.nFaces(), 1.0));
188  scalarField& weight = tweight.ref();
189 
190  // Internal faces
191  forAll(nei, facei)
192  {
193  const point& fc = fCtrs[facei];
194  const vector& fa = fAreas[facei];
195 
196  scalar dOwn = mag(fa & (fc-cellCtrs[own[facei]]));
197  scalar dNei = mag(fa & (cellCtrs[nei[facei]]-fc));
198 
199  weight[facei] = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
200  }
201 
202 
203  // Coupled faces
204 
205  pointField neiCc;
206  syncTools::swapBoundaryCellPositions(mesh, cellCtrs, neiCc);
207 
208  forAll(pbm, patchi)
209  {
210  const polyPatch& pp = pbm[patchi];
211  if (pp.coupled())
212  {
213  forAll(pp, i)
214  {
215  label facei = pp.start() + i;
216  label bFacei = facei - mesh.nInternalFaces();
217 
218  const point& fc = fCtrs[facei];
219  const vector& fa = fAreas[facei];
220 
221  scalar dOwn = mag(fa & (fc-cellCtrs[own[facei]]));
222  scalar dNei = mag(fa & (neiCc[bFacei]-fc));
223 
224  weight[facei] = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
225  }
226  }
227  }
228 
229  return tweight;
230 }
231 
232 
234 (
235  const polyMesh& mesh,
236  const scalarField& vol
237 )
238 {
239  const labelList& own = mesh.faceOwner();
240  const labelList& nei = mesh.faceNeighbour();
241  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
242 
243  tmp<scalarField> tratio(new scalarField(mesh.nFaces(), 1.0));
244  scalarField& ratio = tratio.ref();
245 
246  // Internal faces
247  forAll(nei, facei)
248  {
249  scalar volOwn = vol[own[facei]];
250  scalar volNei = vol[nei[facei]];
251 
252  ratio[facei] = min(volOwn,volNei)/(max(volOwn, volNei)+VSMALL);
253  }
254 
255 
256  // Coupled faces
257 
258  scalarField neiVol;
259  syncTools::swapBoundaryCellList(mesh, vol, neiVol);
260 
261  forAll(pbm, patchi)
262  {
263  const polyPatch& pp = pbm[patchi];
264  if (pp.coupled())
265  {
266  forAll(pp, i)
267  {
268  label facei = pp.start() + i;
269  label bFacei = facei - mesh.nInternalFaces();
270 
271  scalar volOwn = vol[own[facei]];
272  scalar volNei = neiVol[bFacei];
273 
274  ratio[facei] = min(volOwn,volNei)/(max(volOwn, volNei)+VSMALL);
275  }
276  }
277  }
278 
279  return tratio;
280 }
281 
282 
283 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:138
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
primitiveMeshTools.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:329
Foam::polyMeshTools::volRatio
static tmp< scalarField > volRatio(const polyMesh &mesh, const scalarField &vol)
Generate volume ratio field.
Definition: polyMeshTools.C:234
Foam::polyMeshTools::faceSkewness
static tmp< scalarField > faceSkewness(const polyMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
Definition: polyMeshTools.C:91
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::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
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::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:313
polyMeshTools.H
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::polyMeshTools::faceWeights
static tmp< scalarField > faceWeights(const polyMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate interpolation factors field.
Definition: polyMeshTools.C:176
pyramidPointFaceRef.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:36