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-2025 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 
38 {
39  "all",
40  "closest"
41 };
42 
43 namespace Foam
44 {
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::planeToFaceZone::combine(faceZoneSet& fzSet, const bool add) const
53 {
54  // Mark all cells with centres above the plane
55  boolList cellIsAbovePlane(mesh_.nCells());
56  forAll(mesh_.cells(), celli)
57  {
58  cellIsAbovePlane[celli] =
59  ((mesh_.cellCentres()[celli] - point_) & normal_) > 0;
60  }
61 
62  // Mark all coupled neighbour cells with centres above the plane
63  boolList bFaceNbrCellIsAbovePlane(mesh_.nFaces() - mesh_.nInternalFaces());
64  {
65  vectorField bFaceNbrCellCentres;
67  (
68  mesh_,
70  bFaceNbrCellCentres
71  );
72  forAll(bFaceNbrCellIsAbovePlane, bFacei)
73  {
74  bFaceNbrCellIsAbovePlane[bFacei] =
75  ((bFaceNbrCellCentres[bFacei] - point_) & normal_) > 0;
76  }
77  }
78 
79  // Mark all faces that sit between cells above and below the plane
80  boolList faceIsOnPlane(mesh_.nFaces(), false);
81  forAll(mesh_.faceNeighbour(), facei)
82  {
83  faceIsOnPlane[facei] =
84  cellIsAbovePlane[mesh_.faceOwner()[facei]]
85  != cellIsAbovePlane[mesh_.faceNeighbour()[facei]];
86  }
88  {
89  const polyPatch& patch = mesh_.boundaryMesh()[patchi];
90  if (!patch.coupled()) continue;
91  forAll(patch, patchFacei)
92  {
93  const label facei = patch.start() + patchFacei;
94  faceIsOnPlane[facei] =
95  cellIsAbovePlane[mesh_.faceOwner()[facei]]
96  != bFaceNbrCellIsAbovePlane[facei - mesh_.nInternalFaces()];
97  }
98  }
99 
100  // Ensure consistency across couplings
101  syncTools::syncFaceList(mesh_, faceIsOnPlane, orEqOp<bool>());
102 
103  // Convert marked faces to a list of indices
104  labelList newSetFaces(findIndices(faceIsOnPlane, true));
105 
106  // If constructing a single contiguous set, remove all faces except those
107  // connected to the contiguous region closest to the specified point
108  if (include_ == include::closest)
109  {
110  // Step 1: Get locally contiguous regions for the new face set and the
111  // total number of regions across all processors.
112  labelList newSetFaceRegions(newSetFaces.size(), -1);
113  label nRegions = -1;
114  {
115  // Create a patch of the set faces
116  const uindirectPrimitivePatch newSetPatch
117  (
118  UIndirectList<face>(mesh_.faces(), newSetFaces),
119  mesh_.points()
120  );
121 
122  // Get the region ID-s and store the total number of regions on
123  // each processor
124  labelList procNRegions(Pstream::nProcs(), -1);
125  procNRegions[Pstream::myProcNo()] =
127  (
128  newSetPatch,
129  boolList(newSetPatch.nEdges(), false),
130  newSetFaceRegions
131  );
132  Pstream::gatherList(procNRegions);
133  Pstream::scatterList(procNRegions);
134 
135  // Cumulative sum the number of regions on each processor to get an
136  // offset which makes the local region ID-s globally unique
137  labelList procRegionOffset(Pstream::nProcs(), 0);
138  for (label proci = 1; proci < Pstream::nProcs(); ++ proci)
139  {
140  procRegionOffset[proci] +=
141  procRegionOffset[proci - 1]
142  + procNRegions[proci - 1];
143  }
144 
145  // Apply the offset
146  forAll(newSetFaces, newSetFacei)
147  {
148  newSetFaceRegions[newSetFacei] +=
149  procRegionOffset[Pstream::myProcNo()];
150  }
151 
152  // Store the total number of regions across all processors
153  nRegions = procRegionOffset.last() + procNRegions.last();
154  }
155 
156  // Step 2: Create a region map which combines regions which are
157  // connected across coupled interfaces
158  labelList regionMap(identityMap(nRegions));
159  {
160  // Put region labels on connected boundary edges and synchronise to
161  // create a list of all regions connected to a given edge
162  labelListList meshEdgeRegions(mesh_.nEdges(), labelList());
163  forAll(newSetFaces, newSetFacei)
164  {
165  const label facei = newSetFaces[newSetFacei];
166  const label regioni = newSetFaceRegions[newSetFacei];
167  forAll(mesh_.faceEdges()[facei], faceEdgei)
168  {
169  const label edgei = mesh_.faceEdges()[facei][faceEdgei];
170  meshEdgeRegions[edgei] = labelList(1, regioni);
171  }
172  }
174  (
175  mesh_,
176  meshEdgeRegions,
177  globalMeshData::ListPlusEqOp<labelList>(),
178  labelList()
179  );
180 
181  // Combine edge regions to create a list of what regions a given
182  // region is connected to
183  List<labelHashSet> regionRegions(nRegions);
184  forAll(newSetFaces, newSetFacei)
185  {
186  const label facei = newSetFaces[newSetFacei];
187  const label regioni = newSetFaceRegions[newSetFacei];
188  forAll(mesh_.faceEdges()[facei], faceEdgei)
189  {
190  const label edgei = mesh_.faceEdges()[facei][faceEdgei];
191  forAll(meshEdgeRegions[edgei], edgeRegioni)
192  {
193  if (meshEdgeRegions[edgei][edgeRegioni] != regioni)
194  {
195  regionRegions[regioni].insert
196  (
197  meshEdgeRegions[edgei][edgeRegioni]
198  );
199  }
200  }
201  }
202  }
203  Pstream::listCombineGather(regionRegions, plusEqOp<labelHashSet>());
204  Pstream::listCombineScatter(regionRegions);
205 
206  // Collapse the region connections into a map between each region
207  // and the lowest numbered region that it connects to
208  forAll(regionRegions, regioni)
209  {
210  forAllConstIter(labelHashSet, regionRegions[regioni], iter)
211  {
212  regionMap[iter.key()] =
213  min(regionMap[iter.key()], regionMap[regioni]);
214  }
215  }
216  }
217 
218  // Step 3: Combine connected regions
219  labelList regionNFaces;
220  {
221  // Remove duplicates from the region map
222  label regioni0 = 0;
223  forAll(regionMap, regioni)
224  {
225  if (regionMap[regioni] > regioni0)
226  {
227  ++ regioni0;
228  regionMap[regioni] = regioni0;
229  }
230  }
231 
232  // Recompute the number of regions
233  nRegions = regioni0 + 1;
234 
235  // Renumber the face region ID-s
236  newSetFaceRegions =
237  IndirectList<label>(regionMap, newSetFaceRegions);
238 
239  // Report the final number and size of the regions
240  regionNFaces = labelList(nRegions, 0);
241  forAll(newSetFaces, newSetFacei)
242  {
243  regionNFaces[newSetFaceRegions[newSetFacei]] ++;
244  }
245  Pstream::listCombineGather(regionNFaces, plusEqOp<label>());
246  Pstream::listCombineScatter(regionNFaces);
247  Info<< " Found " << nRegions << " contiguous regions with "
248  << regionNFaces << " faces" << endl;
249  }
250 
251  // Step 4: Choose the closest region to output
252  label selectedRegioni = -1;
253  {
254  // Compute the region centres
255  scalarField regionMagAreas(nRegions, 0);
256  pointField regionCentres(nRegions, Zero);
257  forAll(newSetFaces, newSetFacei)
258  {
259  const label facei = newSetFaces[newSetFacei];
260  const label regioni = newSetFaceRegions[newSetFacei];
261 
262  const vector& a = mesh_.faceAreas()[facei];
263  const point& c = mesh_.faceCentres()[facei];
264 
265  regionMagAreas[regioni] += mag(a);
266  regionCentres[regioni] += mag(a)*c;
267  }
268  Pstream::listCombineGather(regionMagAreas, plusEqOp<scalar>());
269  Pstream::listCombineGather(regionCentres, plusEqOp<point>());
270  Pstream::listCombineScatter(regionMagAreas);
271  Pstream::listCombineScatter(regionCentres);
272  regionCentres /= regionMagAreas;
273 
274  // Find the region centroid closest to the reference point
275  selectedRegioni =
277  (
278  findMin(mag(regionCentres - point_)()),
279  minOp<label>()
280  );
281 
282  // Report the selection
283  Info<< " Selecting region " << selectedRegioni << " with "
284  << regionNFaces[selectedRegioni]
285  << " faces as the closest to point " << point_ << endl;
286  }
287 
288  // Step 5: Remove any faces from the set list not in the selected region
289  {
290  // Remove faces from the list by shuffling up and resizing
291  label newSetFacei0 = 0;
292  forAll(newSetFaces, newSetFacei)
293  {
294  newSetFaces[newSetFacei0] = newSetFaces[newSetFacei];
295 
296  if (newSetFaceRegions[newSetFacei] == selectedRegioni)
297  {
298  ++ newSetFacei0;
299  }
300  }
301  newSetFaces.resize(newSetFacei0);
302  }
303  }
304 
305  // Modify the face zone set
306  DynamicList<label> newAddressing;
307  DynamicList<bool> newFlipMap;
308  if (add)
309  {
310  // Start from copy
311  newAddressing = DynamicList<label>(fzSet.addressing());
312  newFlipMap = DynamicList<bool>(fzSet.flipMap());
313 
314  // Add anything from the new set that is not already in the zone set
315  forAll(newSetFaces, newSetFacei)
316  {
317  const label facei = newSetFaces[newSetFacei];
318 
319  if (!fzSet.found(facei))
320  {
321  newAddressing.append(facei);
322  newFlipMap.append(cellIsAbovePlane[mesh_.faceOwner()[facei]]);
323  }
324  }
325  }
326  else
327  {
328  // Start from empty
329  newAddressing = DynamicList<label>(fzSet.addressing().size());
330  newFlipMap = DynamicList<bool>(fzSet.flipMap().size());
331 
332  // Add everything from the zone set that is not also in the new set
333  labelHashSet newSet(newSetFaces);
334  forAll(fzSet.addressing(), i)
335  {
336  const label facei = fzSet.addressing()[i];
337 
338  if (!newSet.found(facei))
339  {
340  newAddressing.append(facei);
341  newFlipMap.append(cellIsAbovePlane[mesh_.faceOwner()[facei]]);
342  }
343  }
344  }
345  fzSet.addressing().transfer(newAddressing);
346  fzSet.flipMap().transfer(newFlipMap);
347  fzSet.updateSet();
348 }
349 
350 
351 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
352 
354 (
355  const polyMesh& mesh,
356  const dictionary& dict
357 )
358 :
360  point_(dict.lookup<vector>("point", dimLength)),
361  normal_(dict.lookup<vector>("normal", dimless)),
362  include_
363  (
364  includeNames_
365  [
366  dict.lookupOrDefault<word>
367  (
368  "include",
369  includeNames_[include::all]
370  )
371  ]
372  )
373 {}
374 
375 
376 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
377 
379 {}
380 
381 
382 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
383 
385 (
386  const topoSetSource::setAction action,
387  topoSet& set
388 ) const
389 {
390  if (!isA<faceZoneSet>(set))
391  {
393  << "Operation only allowed on a faceZoneSet." << endl;
394  }
395  else
396  {
397  faceZoneSet& fzSet = refCast<faceZoneSet>(set);
398 
399  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
400  {
401  Info<< " Adding faces which form a plane at " << point_
402  << " with normal " << normal_ << endl;
403 
404  combine(fzSet, true);
405  }
406  else if (action == topoSetSource::DELETE)
407  {
408  Info<< " Removing faces which form a plane at " << point_
409  << " with normal " << normal_ << endl;
410 
411  combine(fzSet, false);
412  }
413  }
414 }
415 
416 
417 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1363
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
const vectorField & faceCentres() const
label nInternalFaces() const
const vectorField & cellCentres() const
label nCells() const
const labelListList & faceEdges() const
const vectorField & faceAreas() const
const cellList & cells() const
label nEdges() const
label nFaces() 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:61
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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:258
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
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
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)
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
void add(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
treeBoundBox combine(const treeBoundBox &a, const treeBoundBox &b)
Definition: patchToPatch.C:78
dictionary dict