thresholdCellFaces.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "thresholdCellFaces.H"
27 
28 #include "polyMesh.H"
29 #include "DynamicList.H"
30 
31 #include "emptyPolyPatch.H"
32 #include "processorPolyPatch.H"
33 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 defineTypeNameAndDebug(thresholdCellFaces, 0);
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::thresholdCellFaces::calculate
46 (
47  const scalarField& field,
48  const scalar lowerThreshold,
49  const scalar upperThreshold,
50  const bool triangulate
51 )
52 {
53  const labelList& own = mesh_.faceOwner();
54  const labelList& nei = mesh_.faceNeighbour();
55 
56  const faceList& origFaces = mesh_.faces();
57  const pointField& origPoints = mesh_.points();
58 
59  const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
60 
61 
62  surfZoneList surfZones(bMesh.size()+1);
63 
64  surfZones[0] = surfZone
65  (
66  "internalMesh",
67  0, // size
68  0, // start
69  0 // index
70  );
71 
72  forAll(bMesh, patchi)
73  {
74  surfZones[patchi+1] = surfZone
75  (
76  bMesh[patchi].name(),
77  0, // size
78  0, // start
79  patchi+1 // index
80  );
81  }
82 
83 
84  label nFaces = 0;
85  label nPoints = 0;
86 
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  forAll(f, fp)
133  {
134  if (oldToNewPoints[f[fp]] == -1)
135  {
136  oldToNewPoints[f[fp]] = 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  forAll(f, fp)
205  {
206  if (oldToNewPoints[f[fp]] == -1)
207  {
208  oldToNewPoints[f[fp]] = 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  forAll(surfFaces, facei)
241  {
242  inplaceRenumber(oldToNewPoints, surfFaces[facei]);
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 
280  if (lowerThreshold > upperThreshold)
281  {
283  << lowerThreshold << " > " << upperThreshold << endl;
284  }
285 
286  calculate(field, lowerThreshold, upperThreshold, triangulate);
287 }
288 
289 
290 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:59
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
thresholdCellFaces(const polyMesh &, const scalarField &, const scalar lowerThreshold, const scalar upperThreshold, const bool triangulate=false)
Construct from mesh, field and threshold value.
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< label > labelList
A List of labels.
Definition: labelList.H:56
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< surfZone > surfZoneList
Definition: surfZoneList.H:45
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:397
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
label cellId
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Namespace for OpenFOAM.