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-2019 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 #include "emptyPolyPatch.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(thresholdCellFaces, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 void Foam::thresholdCellFaces::calculate
40 (
41  const scalarField& field,
42  const scalar lowerThreshold,
43  const scalar upperThreshold,
44  const bool triangulate
45 )
46 {
47  const labelList& own = mesh_.faceOwner();
48  const labelList& nei = mesh_.faceNeighbour();
49 
50  const faceList& origFaces = mesh_.faces();
51  const pointField& origPoints = mesh_.points();
52 
53  const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
54 
55 
56  surfZoneList surfZones(bMesh.size()+1);
57 
58  surfZones[0] = surfZone
59  (
60  "internalMesh",
61  0, // size
62  0, // start
63  0 // index
64  );
65 
66  forAll(bMesh, patchi)
67  {
68  surfZones[patchi+1] = surfZone
69  (
70  bMesh[patchi].name(),
71  0, // size
72  0, // start
73  patchi+1 // index
74  );
75  }
76 
77 
78  label nFaces = 0;
79  label nPoints = 0;
80 
81 
82  meshCells_.clear();
83 
84  DynamicList<face> surfFaces(0.5 * mesh_.nFaces());
85  DynamicList<label> surfCells(surfFaces.size());
86 
87  labelList oldToNewPoints(origPoints.size(), -1);
88 
89 
90  // internal faces only
91  for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
92  {
93  int side = 0;
94 
95  // check lowerThreshold
96  if (field[own[facei]] > lowerThreshold)
97  {
98  if (field[nei[facei]] < lowerThreshold)
99  {
100  side = +1;
101  }
102  }
103  else if (field[nei[facei]] > lowerThreshold)
104  {
105  side = -1;
106  }
107 
108  // check upperThreshold
109  if (field[own[facei]] < upperThreshold)
110  {
111  if (field[nei[facei]] > upperThreshold)
112  {
113  side = +1;
114  }
115  }
116  else if (field[nei[facei]] < upperThreshold)
117  {
118  side = -1;
119  }
120 
121 
122  if (side)
123  {
124  const face& f = origFaces[facei];
125 
126  forAll(f, fp)
127  {
128  if (oldToNewPoints[f[fp]] == -1)
129  {
130  oldToNewPoints[f[fp]] = nPoints++;
131  }
132  }
133 
134 
135  label cellId;
136  face surfFace;
137 
138  if (side > 0)
139  {
140  surfFace = f;
141  cellId = own[facei];
142  }
143  else
144  {
145  surfFace = f.reverseFace();
146  cellId = nei[facei];
147  }
148 
149 
150  if (triangulate)
151  {
152  label count = surfFace.triangles(origPoints, surfFaces);
153  while (count-- > 0)
154  {
155  surfCells.append(cellId);
156  }
157  }
158  else
159  {
160  surfFaces.append(surfFace);
161  surfCells.append(cellId);
162  }
163  }
164  }
165 
166  surfZones[0].size() = surfFaces.size();
167 
168 
169  // nothing special for processor patches?
170  forAll(bMesh, patchi)
171  {
172  const polyPatch& p = bMesh[patchi];
173  surfZone& zone = surfZones[patchi+1];
174 
175  zone.start() = nFaces;
176 
177  if
178  (
179  isA<emptyPolyPatch>(p)
180  || (Pstream::parRun() && isA<processorPolyPatch>(p))
181  )
182  {
183  continue;
184  }
185 
186  label facei = p.start();
187 
188  // patch faces
189  forAll(p, localFacei)
190  {
191  if
192  (
193  field[own[facei]] > lowerThreshold
194  && field[own[facei]] < upperThreshold
195  )
196  {
197  const face& f = origFaces[facei];
198  forAll(f, fp)
199  {
200  if (oldToNewPoints[f[fp]] == -1)
201  {
202  oldToNewPoints[f[fp]] = nPoints++;
203  }
204  }
205 
206  label cellId = own[facei];
207 
208  if (triangulate)
209  {
210  label count = f.triangles(origPoints, surfFaces);
211  while (count-- > 0)
212  {
213  surfCells.append(cellId);
214  }
215  }
216  else
217  {
218  surfFaces.append(f);
219  surfCells.append(cellId);
220  }
221  }
222 
223  ++facei;
224  }
225 
226  zone.size() = surfFaces.size() - zone.start();
227  }
228 
229 
230  surfFaces.shrink();
231  surfCells.shrink();
232 
233  // renumber
234  forAll(surfFaces, facei)
235  {
236  inplaceRenumber(oldToNewPoints, surfFaces[facei]);
237  }
238 
239 
240  pointField surfPoints(nPoints);
241  nPoints = 0;
242  forAll(oldToNewPoints, pointi)
243  {
244  if (oldToNewPoints[pointi] >= 0)
245  {
246  surfPoints[oldToNewPoints[pointi]] = origPoints[pointi];
247  nPoints++;
248  }
249  }
250  surfPoints.setSize(nPoints);
251 
252  this->storedPoints().transfer(surfPoints);
253  this->storedFaces().transfer(surfFaces);
254  this->storedZones().transfer(surfZones);
255 
256  meshCells_.transfer(surfCells);
257 }
258 
259 
260 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
261 
263 (
264  const polyMesh& mesh,
265  const scalarField& field,
266  const scalar lowerThreshold,
267  const scalar upperThreshold,
268  const bool triangulate
269 )
270 :
271  mesh_(mesh)
272 {
273 
274  if (lowerThreshold > upperThreshold)
275  {
277  << lowerThreshold << " > " << upperThreshold << endl;
278  }
279 
280  calculate(field, lowerThreshold, upperThreshold, triangulate);
281 }
282 
283 
284 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:251
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:399
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.