Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
34
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
41}
42
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
46void 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
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 (
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);
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// ************************************************************************* //
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
const surfZoneList & surfZones() const
virtual label triangulate()
const List< face > & surfFaces() const
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
Selects the mesh cell faces specified by a threshold value. Non-triangulated by default.
thresholdCellFaces(const polyMesh &mesh, const scalarField &field, const scalar lowerThreshold, const scalar upperThreshold, const bool triangulate=false)
Construct from mesh, field and threshold values.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
rDeltaTY field()
dynamicFvMesh & mesh
auto & name
label cellId
#define WarningInFunction
Report a warning using Foam::Warning.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition BitOps.H:73
Namespace for OpenFOAM.
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
List< label > labelList
A List of labels.
Definition List.H:62
List< surfZone > surfZoneList
List of surfZone.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points).
Definition bMesh.H:39
vectorField pointField
pointField is a vectorField.
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299