thresholdCellFaces.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "thresholdCellFaces.H"
30 
31 #include "polyMesh.H"
32 #include "DynamicList.H"
33 
34 #include "emptyPolyPatch.H"
35 #include "processorPolyPatch.H"
36 
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 defineTypeNameAndDebug(thresholdCellFaces, 0);
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::thresholdCellFaces::calculate
49 (
50  const scalarField& field,
51  const scalar lowerThreshold,
52  const scalar upperThreshold,
53  const bool triangulate
54 )
55 {
56  const labelList& own = mesh_.faceOwner();
57  const labelList& nei = mesh_.faceNeighbour();
58 
59  const faceList& origFaces = mesh_.faces();
60  const pointField& origPoints = mesh_.points();
61 
62  const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
63 
64 
65  surfZoneList surfZones(bMesh.size()+1);
66 
67  surfZones[0] = surfZone
68  (
69  "internalMesh",
70  0, // size
71  0, // start
72  0 // index
73  );
74 
75  forAll(bMesh, patchi)
76  {
77  surfZones[patchi+1] = surfZone
78  (
79  bMesh[patchi].name(),
80  0, // size
81  0, // start
82  patchi+1 // index
83  );
84  }
85 
86 
87  label nFaces = 0;
88  label nPoints = 0;
89 
90  meshCells_.clear();
91 
92  DynamicList<face> surfFaces(0.5 * mesh_.nFaces());
93  DynamicList<label> surfCells(surfFaces.size());
94 
95  labelList oldToNewPoints(origPoints.size(), -1);
96 
97 
98  // internal faces only
99  for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
100  {
101  int side = 0;
102 
103  // Check lowerThreshold
104  if (field[own[facei]] > lowerThreshold)
105  {
106  if (field[nei[facei]] < lowerThreshold)
107  {
108  side = +1;
109  }
110  }
111  else if (field[nei[facei]] > lowerThreshold)
112  {
113  side = -1;
114  }
115 
116  // Check upperThreshold
117  if (field[own[facei]] < upperThreshold)
118  {
119  if (field[nei[facei]] > upperThreshold)
120  {
121  side = +1;
122  }
123  }
124  else if (field[nei[facei]] < upperThreshold)
125  {
126  side = -1;
127  }
128 
129 
130  if (side)
131  {
132  const face& f = origFaces[facei];
133 
134  for (const label pointi : f)
135  {
136  if (oldToNewPoints[pointi] == -1)
137  {
138  oldToNewPoints[pointi] = nPoints++;
139  }
140  }
141 
142 
143  label cellId;
144  face surfFace;
145 
146  if (side > 0)
147  {
148  surfFace = f;
149  cellId = own[facei];
150  }
151  else
152  {
153  surfFace = f.reverseFace();
154  cellId = nei[facei];
155  }
156 
157 
158  if (triangulate)
159  {
160  label count = surfFace.triangles(origPoints, surfFaces);
161  while (count-- > 0)
162  {
163  surfCells.append(cellId);
164  }
165  }
166  else
167  {
168  surfFaces.append(surfFace);
169  surfCells.append(cellId);
170  }
171  }
172  }
173 
174  surfZones[0].size() = surfFaces.size();
175 
176 
177  // nothing special for processor patches?
178  forAll(bMesh, patchi)
179  {
180  const polyPatch& p = bMesh[patchi];
181  surfZone& zone = surfZones[patchi+1];
182 
183  zone.start() = nFaces;
184 
185  if
186  (
187  isA<emptyPolyPatch>(p)
188  || (Pstream::parRun() && isA<processorPolyPatch>(p))
189  )
190  {
191  continue;
192  }
193 
194  label facei = p.start();
195 
196  // patch faces
197  forAll(p, localFacei)
198  {
199  if
200  (
201  field[own[facei]] > lowerThreshold
202  && field[own[facei]] < upperThreshold
203  )
204  {
205  const face& f = origFaces[facei];
206  for (const label pointi : f)
207  {
208  if (oldToNewPoints[pointi] == -1)
209  {
210  oldToNewPoints[pointi] = nPoints++;
211  }
212  }
213 
214  label cellId = own[facei];
215 
216  if (triangulate)
217  {
218  label count = f.triangles(origPoints, surfFaces);
219  while (count-- > 0)
220  {
221  surfCells.append(cellId);
222  }
223  }
224  else
225  {
226  surfFaces.append(f);
227  surfCells.append(cellId);
228  }
229  }
230 
231  ++facei;
232  }
233 
234  zone.size() = surfFaces.size() - zone.start();
235  }
236 
237 
238  surfFaces.shrink();
239  surfCells.shrink();
240 
241  // renumber
242  for (face& f : surfFaces)
243  {
244  inplaceRenumber(oldToNewPoints, f);
245  }
246 
247 
248  pointField surfPoints(nPoints);
249  nPoints = 0;
250  forAll(oldToNewPoints, pointi)
251  {
252  if (oldToNewPoints[pointi] >= 0)
253  {
254  surfPoints[oldToNewPoints[pointi]] = origPoints[pointi];
255  ++nPoints;
256  }
257  }
258  surfPoints.setSize(nPoints);
259 
260  this->storedPoints().transfer(surfPoints);
261  this->storedFaces().transfer(surfFaces);
262  this->storedZones().transfer(surfZones);
263 
264  meshCells_.transfer(surfCells);
265 }
266 
267 
268 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
269 
271 (
272  const polyMesh& mesh,
273  const scalarField& field,
274  const scalar lowerThreshold,
275  const scalar upperThreshold,
276  const bool triangulate
277 )
278 :
279  mesh_(mesh)
280 {
281  if (lowerThreshold > upperThreshold)
282  {
284  << lowerThreshold << " > " << upperThreshold << endl;
285  }
286 
287  calculate(field, lowerThreshold, upperThreshold, triangulate);
288 }
289 
290 
291 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
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:1038
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::MeshedSurface< face >::surfZones
const surfZoneList & surfZones() const
Const access to the surface zones.
Definition: MeshedSurface.H:370
Foam::MeshedSurface< face >::storedFaces
List< face > & storedFaces()
Non-const access to the faces.
Definition: MeshedSurface.H:157
Foam::MeshedSurface< face >::triangulate
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
Definition: MeshedSurface.C:864
Foam::thresholdCellFaces::thresholdCellFaces
thresholdCellFaces(const polyMesh &mesh, const scalarField &field, const scalar lowerThreshold, const scalar upperThreshold, const bool triangulate=false)
Construct from mesh, field and threshold values.
Definition: thresholdCellFaces.C:271
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::MeshedSurface< face >::storedPoints
pointField & storedPoints()
Non-const access to global points.
Definition: MeshedSurface.H:151
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::MeshedSurface< face >::surfFaces
const List< face > & surfFaces() const
Return const access to the faces.
Definition: MeshedSurface.H:362
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
polyMesh.H
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:61
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:290
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::MeshedSurface< face >::storedZones
surfZoneList & storedZones()
Non-const access to the zones.
Definition: MeshedSurface.H:163
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< scalar >
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::bMesh
PrimitivePatch< face, List, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points)
Definition: bMesh.H:47
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1076
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:436
field
rDeltaTY field()
thresholdCellFaces.H
processorPolyPatch.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
emptyPolyPatch.H
cellId
label cellId
Definition: interrogateWallPatches.H:67
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
f
labelList f(nPoints)
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:74
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
Foam::surfZoneList
List< surfZone > surfZoneList
Definition: surfZoneList.H:47
DynamicList.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1082