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 #include "polyMesh.H"
31 #include "DynamicList.H"
32 #include "emptyPolyPatch.H"
33 #include "processorPolyPatch.H"
34 
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 defineTypeNameAndDebug(thresholdCellFaces, 0);
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::thresholdCellFaces::calculate
47 (
48  const scalarField& field,
49  const scalar lowerThreshold,
50  const scalar upperThreshold,
51  const bool triangulate
52 )
53 {
54  const labelList& own = mesh_.faceOwner();
55  const labelList& nei = mesh_.faceNeighbour();
56 
57  const faceList& origFaces = mesh_.faces();
58  const pointField& origPoints = mesh_.points();
59 
60  const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
61 
62 
63  surfZoneList surfZones(bMesh.size()+1);
64 
65  surfZones[0] = surfZone
66  (
67  "internalMesh",
68  0, // size
69  0, // start
70  0 // index
71  );
72 
73  forAll(bMesh, patchi)
74  {
75  surfZones[patchi+1] = surfZone
76  (
77  bMesh[patchi].name(),
78  0, // size
79  0, // start
80  patchi+1 // index
81  );
82  }
83 
84 
85  label nFaces = 0;
86  label nPoints = 0;
87 
88  meshCells_.clear();
89 
90  DynamicList<face> surfFaces(0.5 * mesh_.nFaces());
91  DynamicList<label> surfCells(surfFaces.size());
92 
93  labelList oldToNewPoints(origPoints.size(), -1);
94 
95 
96  // internal faces only
97  for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
98  {
99  int side = 0;
100 
101  // Check lowerThreshold
102  if (field[own[facei]] > lowerThreshold)
103  {
104  if (field[nei[facei]] < lowerThreshold)
105  {
106  side = +1;
107  }
108  }
109  else if (field[nei[facei]] > lowerThreshold)
110  {
111  side = -1;
112  }
113 
114  // Check upperThreshold
115  if (field[own[facei]] < upperThreshold)
116  {
117  if (field[nei[facei]] > upperThreshold)
118  {
119  side = +1;
120  }
121  }
122  else if (field[nei[facei]] < upperThreshold)
123  {
124  side = -1;
125  }
126 
127 
128  if (side)
129  {
130  const face& f = origFaces[facei];
131 
132  for (const label pointi : f)
133  {
134  if (oldToNewPoints[pointi] == -1)
135  {
136  oldToNewPoints[pointi] = nPoints++;
137  }
138  }
139 
140 
141  label cellId;
142  face surfFace;
143 
144  if (side > 0)
145  {
146  surfFace = f;
147  cellId = own[facei];
148  }
149  else
150  {
151  surfFace = f.reverseFace();
152  cellId = nei[facei];
153  }
154 
155 
156  if (triangulate)
157  {
158  label count = surfFace.triangles(origPoints, surfFaces);
159  while (count-- > 0)
160  {
161  surfCells.append(cellId);
162  }
163  }
164  else
165  {
166  surfFaces.append(surfFace);
167  surfCells.append(cellId);
168  }
169  }
170  }
171 
172  surfZones[0].size() = surfFaces.size();
173 
174 
175  // nothing special for processor patches?
176  forAll(bMesh, patchi)
177  {
178  const polyPatch& p = bMesh[patchi];
179  surfZone& zone = surfZones[patchi+1];
180 
181  zone.start() = nFaces;
182 
183  if
184  (
185  isA<emptyPolyPatch>(p)
186  || (Pstream::parRun() && isA<processorPolyPatch>(p))
187  )
188  {
189  continue;
190  }
191 
192  label facei = p.start();
193 
194  // patch faces
195  forAll(p, localFacei)
196  {
197  if
198  (
199  field[own[facei]] > lowerThreshold
200  && field[own[facei]] < upperThreshold
201  )
202  {
203  const face& f = origFaces[facei];
204  for (const label pointi : f)
205  {
206  if (oldToNewPoints[pointi] == -1)
207  {
208  oldToNewPoints[pointi] = nPoints++;
209  }
210  }
211 
212  label cellId = own[facei];
213 
214  if (triangulate)
215  {
216  label count = f.triangles(origPoints, surfFaces);
217  while (count-- > 0)
218  {
219  surfCells.append(cellId);
220  }
221  }
222  else
223  {
224  surfFaces.append(f);
225  surfCells.append(cellId);
226  }
227  }
228 
229  ++facei;
230  }
231 
232  zone.size() = surfFaces.size() - zone.start();
233  }
234 
235 
236  surfFaces.shrink();
237  surfCells.shrink();
238 
239  // renumber
240  for (face& f : surfFaces)
241  {
242  inplaceRenumber(oldToNewPoints, f);
243  }
244 
245 
246  pointField surfPoints(nPoints);
247  nPoints = 0;
248  forAll(oldToNewPoints, pointi)
249  {
250  if (oldToNewPoints[pointi] >= 0)
251  {
252  surfPoints[oldToNewPoints[pointi]] = origPoints[pointi];
253  ++nPoints;
254  }
255  }
256  surfPoints.setSize(nPoints);
257 
258  this->storedPoints().transfer(surfPoints);
259  this->storedFaces().transfer(surfFaces);
260  this->storedZones().transfer(surfZones);
261 
262  meshCells_.transfer(surfCells);
263 }
264 
265 
266 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
267 
269 (
270  const polyMesh& mesh,
271  const scalarField& field,
272  const scalar lowerThreshold,
273  const scalar upperThreshold,
274  const bool triangulate
275 )
276 :
277  mesh_(mesh)
278 {
279  if (lowerThreshold > upperThreshold)
280  {
282  << lowerThreshold << " > " << upperThreshold << endl;
283  }
284 
285  calculate(field, lowerThreshold, upperThreshold, triangulate);
286 }
287 
288 
289 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
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
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::MeshedSurface::surfZones
const surfZoneList & surfZones() const
Const access to the surface zones.
Definition: MeshedSurface.H:429
Foam::MeshedSurface::storedFaces
List< Face > & storedFaces()
Non-const access to the faces.
Definition: MeshedSurface.H:173
Foam::MeshedSurface::triangulate
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
Definition: MeshedSurface.C:1015
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:269
Foam::MeshedSurface::storedPoints
pointField & storedPoints()
Non-const access to global points.
Definition: MeshedSurface.H:167
Foam::MeshedSurface::surfFaces
const List< Face > & surfFaces() const
Return const access to the faces.
Definition: MeshedSurface.H:413
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
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:296
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::MeshedSurface::storedZones
surfZoneList & storedZones()
Non-const access to the zones.
Definition: MeshedSurface.H:185
Foam::Field< scalar >
Foam::bMesh
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points)
Definition: bMesh.H:48
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
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:1094
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:77
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::surfZoneList
List< surfZone > surfZoneList
Definition: surfZoneList.H:47
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
DynamicList.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113