planeToFaceZone.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) 2020-2024 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 "planeToFaceZone.H"
27 #include "polyMesh.H"
28 #include "faceZoneSet.H"
30 #include "PatchTools.H"
31 #include "syncTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  template<>
40  {
41  "all",
42  "closest"
43  };
44 }
45 
48 
49 
50 namespace Foam
51 {
54 }
55 
56 
57 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58 
59 void Foam::planeToFaceZone::combine(faceZoneSet& fzSet, const bool add) const
60 {
61  // Mark all cells with centres above the plane
62  boolList cellIsAbovePlane(mesh_.nCells());
63  forAll(mesh_.cells(), celli)
64  {
65  cellIsAbovePlane[celli] =
66  ((mesh_.cellCentres()[celli] - point_) & normal_) > 0;
67  }
68 
69  // Mark all coupled neighbour cells with centres above the plane
70  boolList bFaceNbrCellIsAbovePlane(mesh_.nFaces() - mesh_.nInternalFaces());
71  {
72  vectorField bFaceNbrCellCentres;
74  (
75  mesh_,
77  bFaceNbrCellCentres
78  );
79  forAll(bFaceNbrCellIsAbovePlane, bFacei)
80  {
81  bFaceNbrCellIsAbovePlane[bFacei] =
82  ((bFaceNbrCellCentres[bFacei] - point_) & normal_) > 0;
83  }
84  }
85 
86  // Mark all faces that sit between cells above and below the plane
87  boolList faceIsOnPlane(mesh_.nFaces(), false);
88  forAll(mesh_.faceNeighbour(), facei)
89  {
90  faceIsOnPlane[facei] =
91  cellIsAbovePlane[mesh_.faceOwner()[facei]]
92  != cellIsAbovePlane[mesh_.faceNeighbour()[facei]];
93  }
95  {
96  const polyPatch& patch = mesh_.boundaryMesh()[patchi];
97  if (!patch.coupled()) continue;
98  forAll(patch, patchFacei)
99  {
100  const label facei = patch.start() + patchFacei;
101  faceIsOnPlane[facei] =
102  cellIsAbovePlane[mesh_.faceOwner()[facei]]
103  != bFaceNbrCellIsAbovePlane[facei - mesh_.nInternalFaces()];
104  }
105  }
106 
107  // Ensure consistency across couplings
108  syncTools::syncFaceList(mesh_, faceIsOnPlane, orEqOp<bool>());
109 
110  // Convert marked faces to a list of indices
111  labelList newSetFaces(findIndices(faceIsOnPlane, true));
112 
113  // If constructing a single contiguous set, remove all faces except those
114  // connected to the contiguous region closest to the specified point
115  if (include_ == include::closest)
116  {
117  // Step 1: Get locally contiguous regions for the new face set and the
118  // total number of regions across all processors.
119  labelList newSetFaceRegions(newSetFaces.size(), -1);
120  label nRegions = -1;
121  {
122  // Create a patch of the set faces
123  const uindirectPrimitivePatch newSetPatch
124  (
125  UIndirectList<face>(mesh_.faces(), newSetFaces),
126  mesh_.points()
127  );
128 
129  // Get the region ID-s and store the total number of regions on
130  // each processor
131  labelList procNRegions(Pstream::nProcs(), -1);
132  procNRegions[Pstream::myProcNo()] =
134  (
135  newSetPatch,
136  boolList(newSetPatch.nEdges(), false),
137  newSetFaceRegions
138  );
139  Pstream::gatherList(procNRegions);
140  Pstream::scatterList(procNRegions);
141 
142  // Cumulative sum the number of regions on each processor to get an
143  // offset which makes the local region ID-s globally unique
144  labelList procRegionOffset(Pstream::nProcs(), 0);
145  for (label proci = 1; proci < Pstream::nProcs(); ++ proci)
146  {
147  procRegionOffset[proci] +=
148  procRegionOffset[proci - 1]
149  + procNRegions[proci - 1];
150  }
151 
152  // Apply the offset
153  forAll(newSetFaces, newSetFacei)
154  {
155  newSetFaceRegions[newSetFacei] +=
156  procRegionOffset[Pstream::myProcNo()];
157  }
158 
159  // Store the total number of regions across all processors
160  nRegions = procRegionOffset.last() + procNRegions.last();
161  }
162 
163  // Step 2: Create a region map which combines regions which are
164  // connected across coupled interfaces
165  labelList regionMap(identityMap(nRegions));
166  {
167  // Put region labels on connected boundary edges and synchronise to
168  // create a list of all regions connected to a given edge
169  labelListList meshEdgeRegions(mesh_.nEdges(), labelList());
170  forAll(newSetFaces, newSetFacei)
171  {
172  const label facei = newSetFaces[newSetFacei];
173  const label regioni = newSetFaceRegions[newSetFacei];
174  forAll(mesh_.faceEdges()[facei], faceEdgei)
175  {
176  const label edgei = mesh_.faceEdges()[facei][faceEdgei];
177  meshEdgeRegions[edgei] = labelList(1, regioni);
178  }
179  }
181  (
182  mesh_,
183  meshEdgeRegions,
184  globalMeshData::ListPlusEqOp<labelList>(),
185  labelList()
186  );
187 
188  // Combine edge regions to create a list of what regions a given
189  // region is connected to
190  List<labelHashSet> regionRegions(nRegions);
191  forAll(newSetFaces, newSetFacei)
192  {
193  const label facei = newSetFaces[newSetFacei];
194  const label regioni = newSetFaceRegions[newSetFacei];
195  forAll(mesh_.faceEdges()[facei], faceEdgei)
196  {
197  const label edgei = mesh_.faceEdges()[facei][faceEdgei];
198  forAll(meshEdgeRegions[edgei], edgeRegioni)
199  {
200  if (meshEdgeRegions[edgei][edgeRegioni] != regioni)
201  {
202  regionRegions[regioni].insert
203  (
204  meshEdgeRegions[edgei][edgeRegioni]
205  );
206  }
207  }
208  }
209  }
210  Pstream::listCombineGather(regionRegions, plusEqOp<labelHashSet>());
211  Pstream::listCombineScatter(regionRegions);
212 
213  // Collapse the region connections into a map between each region
214  // and the lowest numbered region that it connects to
215  forAll(regionRegions, regioni)
216  {
217  forAllConstIter(labelHashSet, regionRegions[regioni], iter)
218  {
219  regionMap[iter.key()] =
220  min(regionMap[iter.key()], regionMap[regioni]);
221  }
222  }
223  }
224 
225  // Step 3: Combine connected regions
226  labelList regionNFaces;
227  {
228  // Remove duplicates from the region map
229  label regioni0 = 0;
230  forAll(regionMap, regioni)
231  {
232  if (regionMap[regioni] > regioni0)
233  {
234  ++ regioni0;
235  regionMap[regioni] = regioni0;
236  }
237  }
238 
239  // Recompute the number of regions
240  nRegions = regioni0 + 1;
241 
242  // Renumber the face region ID-s
243  newSetFaceRegions =
244  IndirectList<label>(regionMap, newSetFaceRegions);
245 
246  // Report the final number and size of the regions
247  regionNFaces = labelList(nRegions, 0);
248  forAll(newSetFaces, newSetFacei)
249  {
250  regionNFaces[newSetFaceRegions[newSetFacei]] ++;
251  }
252  Pstream::listCombineGather(regionNFaces, plusEqOp<label>());
253  Pstream::listCombineScatter(regionNFaces);
254  Info<< " Found " << nRegions << " contiguous regions with "
255  << regionNFaces << " faces" << endl;
256  }
257 
258  // Step 4: Choose the closest region to output
259  label selectedRegioni = -1;
260  {
261  // Compute the region centres
262  scalarField regionMagAreas(nRegions, 0);
263  pointField regionCentres(nRegions, Zero);
264  forAll(newSetFaces, newSetFacei)
265  {
266  const label facei = newSetFaces[newSetFacei];
267  const label regioni = newSetFaceRegions[newSetFacei];
268 
269  const vector& a = mesh_.faceAreas()[facei];
270  const point& c = mesh_.faceCentres()[facei];
271 
272  regionMagAreas[regioni] += mag(a);
273  regionCentres[regioni] += mag(a)*c;
274  }
275  Pstream::listCombineGather(regionMagAreas, plusEqOp<scalar>());
276  Pstream::listCombineGather(regionCentres, plusEqOp<point>());
277  Pstream::listCombineScatter(regionMagAreas);
278  Pstream::listCombineScatter(regionCentres);
279  regionCentres /= regionMagAreas;
280 
281  // Find the region centroid closest to the reference point
282  selectedRegioni =
284  (
285  findMin(mag(regionCentres - point_)()),
286  minOp<label>()
287  );
288 
289  // Report the selection
290  Info<< " Selecting region " << selectedRegioni << " with "
291  << regionNFaces[selectedRegioni]
292  << " faces as the closest to point " << point_ << endl;
293  }
294 
295  // Step 5: Remove any faces from the set list not in the selected region
296  {
297  // Remove faces from the list by shuffling up and resizing
298  label newSetFacei0 = 0;
299  forAll(newSetFaces, newSetFacei)
300  {
301  newSetFaces[newSetFacei0] = newSetFaces[newSetFacei];
302 
303  if (newSetFaceRegions[newSetFacei] == selectedRegioni)
304  {
305  ++ newSetFacei0;
306  }
307  }
308  newSetFaces.resize(newSetFacei0);
309  }
310  }
311 
312  // Modify the face zone set
313  DynamicList<label> newAddressing;
314  DynamicList<bool> newFlipMap;
315  if (add)
316  {
317  // Start from copy
318  newAddressing = DynamicList<label>(fzSet.addressing());
319  newFlipMap = DynamicList<bool>(fzSet.flipMap());
320 
321  // Add anything from the new set that is not already in the zone set
322  forAll(newSetFaces, newSetFacei)
323  {
324  const label facei = newSetFaces[newSetFacei];
325 
326  if (!fzSet.found(facei))
327  {
328  newAddressing.append(facei);
329  newFlipMap.append(cellIsAbovePlane[mesh_.faceOwner()[facei]]);
330  }
331  }
332  }
333  else
334  {
335  // Start from empty
336  newAddressing = DynamicList<label>(fzSet.addressing().size());
337  newFlipMap = DynamicList<bool>(fzSet.flipMap().size());
338 
339  // Add everything from the zone set that is not also in the new set
340  labelHashSet newSet(newSetFaces);
341  forAll(fzSet.addressing(), i)
342  {
343  const label facei = fzSet.addressing()[i];
344 
345  if (!newSet.found(facei))
346  {
347  newAddressing.append(facei);
348  newFlipMap.append(cellIsAbovePlane[mesh_.faceOwner()[facei]]);
349  }
350  }
351  }
352  fzSet.addressing().transfer(newAddressing);
353  fzSet.flipMap().transfer(newFlipMap);
354  fzSet.updateSet();
355 }
356 
357 
358 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
359 
361 (
362  const polyMesh& mesh,
363  const dictionary& dict
364 )
365 :
366  topoSetSource(mesh),
367  point_(dict.lookup<vector>("point", dimLength)),
368  normal_(dict.lookup<vector>("normal", dimless)),
369  include_
370  (
371  includeNames_
372  [
373  dict.lookupOrDefault<word>
374  (
375  "include",
376  includeNames_[include::all]
377  )
378  ]
379  )
380 {}
381 
382 
383 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
384 
386 {}
387 
388 
389 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
390 
392 (
393  const topoSetSource::setAction action,
394  topoSet& set
395 ) const
396 {
397  if (!isA<faceZoneSet>(set))
398  {
400  << "Operation only allowed on a faceZoneSet." << endl;
401  }
402  else
403  {
404  faceZoneSet& fzSet = refCast<faceZoneSet>(set);
405 
406  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
407  {
408  Info<< " Adding faces which form a plane at " << point_
409  << " with normal " << normal_ << endl;
410 
411  combine(fzSet, true);
412  }
413  else if (action == topoSetSource::DELETE)
414  {
415  Info<< " Removing faces which form a plane at " << point_
416  << " with normal " << normal_ << endl;
417 
418  combine(fzSet, false);
419  }
420  }
421 }
422 
423 
424 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
static label markZones(const PrimitivePatch< FaceList, PointField > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Like faceSet but -reads data from faceZone -updates faceZone when writing.
Definition: faceZoneSet.H:52
A topoSetSource to select faces based on the adjacent cell centres spanning a given plane....
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
static const NamedEnum< include, 2 > includeNames_
Included region names.
include
Enumeration for what to include.
virtual ~planeToFaceZone()
Destructor.
planeToFaceZone(const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1345
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
label nEdges() const
label nCells() const
const vectorField & faceCentres() const
const vectorField & cellCentres() const
label nInternalFaces() const
label nFaces() const
const labelListList & faceEdges() const
const vectorField & faceAreas() const
const cellList & cells() const
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:387
static void swapBoundaryCellPositions(const polyMesh &mesh, const UList< point > &cellData, List< point > &neighbourCellData)
Swap to obtain neighbour cell positions for all boundary faces.
Definition: syncTools.C:31
Base class of a source for a topoSet.
Definition: topoSetSource.H:64
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:83
const polyMesh & mesh_
Definition: topoSetSource.H:99
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:65
A class for handling words, derived from string.
Definition: word.H:62
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
vector point
Point is a vector.
Definition: point.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
Foam::uindirectPrimitivePatch.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
Field< vector > vectorField
Specialisation of Field<T> for vector.
label findMin(const ListType &, const label start=0)
Find index of min element (and less than given element).
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
treeBoundBox combine(const treeBoundBox &a, const treeBoundBox &b)
Definition: patchToPatch.C:78
dictionary dict