planeToFaceZone.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) 2020 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "planeToFaceZone.H"
30 #include "polyMesh.H"
31 #include "faceZoneSet.H"
33 #include "PatchTools.H"
34 #include "syncTools.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(planeToFaceZone, 0);
42  addToRunTimeSelectionTable(topoSetSource, planeToFaceZone, word);
43  addToRunTimeSelectionTable(topoSetSource, planeToFaceZone, istream);
44 
45  addToRunTimeSelectionTable(topoSetFaceZoneSource, planeToFaceZone, word);
46  addToRunTimeSelectionTable(topoSetFaceZoneSource, planeToFaceZone, istream);
47 
49  (
50  topoSetFaceZoneSource,
51  planeToFaceZone,
52  word,
53  plane
54  );
56  (
57  topoSetFaceZoneSource,
58  planeToFaceZone,
59  istream,
60  plane
61  );
62 }
63 
64 
65 Foam::topoSetSource::addToUsageTable Foam::planeToFaceZone::usage_
66 (
67  planeToFaceZone::typeName,
68  "\n Usage: planeToFaceZone (px py pz) (nx ny nz) include\n\n"
69  " Select faces for which the adjacent cell centres lie on opposite "
70  " of a plane\n\n"
71 );
72 
73 
74 const Foam::Enum
75 <
77 >
78 Foam::planeToFaceZone::faceActionNames_
79 ({
80  { faceAction::ALL, "all" },
81  { faceAction::CLOSEST, "closest" },
82 });
83 
84 
85 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
86 
87 void Foam::planeToFaceZone::combine(faceZoneSet& fzSet, const bool add) const
88 {
89  // Mark all cells with centres above the plane
90  bitSet cellIsAbovePlane(mesh_.nCells());
91  forAll(mesh_.cells(), celli)
92  {
93  if (((mesh_.cellCentres()[celli] - point_) & normal_) > 0)
94  {
95  cellIsAbovePlane.set(celli);
96  }
97  }
98 
99  // Mark all faces that sit between cells above and below the plane
100  bitSet faceIsOnPlane(mesh_.nFaces());
101  forAll(mesh_.faceNeighbour(), facei)
102  {
103  if
104  (
105  cellIsAbovePlane[mesh_.faceOwner()[facei]]
106  != cellIsAbovePlane[mesh_.faceNeighbour()[facei]]
107  )
108  {
109  faceIsOnPlane.set(facei);
110  }
111  }
112  forAll(mesh_.boundaryMesh(), patchi)
113  {
114  const polyPatch& patch = mesh_.boundaryMesh()[patchi];
115  forAll(patch, patchFacei)
116  {
117  const label facei = patch.start() + patchFacei;
118  if (patch.coupled() && cellIsAbovePlane[mesh_.faceOwner()[facei]])
119  {
120  faceIsOnPlane.set(facei);
121  }
122  }
123  }
124  syncTools::syncFaceList(mesh_, faceIsOnPlane, xorEqOp<unsigned int>());
125 
126  // Convert marked faces to a list of indices
127  labelList newSetFaces(faceIsOnPlane.sortedToc());
128  faceIsOnPlane.clear();
129 
130  // If constructing a single contiguous set, remove all faces except those
131  // connected to the contiguous region closest to the specified point
132  if (option_ == faceAction::CLOSEST)
133  {
134  // Step 1: Get locally contiguous regions for the new face set and the
135  // total number of regions across all processors.
136  labelList newSetFaceRegions(newSetFaces.size(), -1);
137  label nRegions = -1;
138  {
139  // Create a patch of the set faces
140  const uindirectPrimitivePatch newSetPatch
141  (
142  UIndirectList<face>(mesh_.faces(), newSetFaces),
143  mesh_.points()
144  );
145 
146  // Get the region ID-s and store the total number of regions on
147  // each processor
148  labelList procNRegions(Pstream::nProcs(), -1);
149  procNRegions[Pstream::myProcNo()] =
151  (
152  newSetPatch,
153  bitSet(), // No border edges
154  newSetFaceRegions
155  );
156  Pstream::gatherList(procNRegions);
157  Pstream::scatterList(procNRegions);
158 
159  // Cumulative sum the number of regions on each processor to get an
160  // offset which makes the local region ID-s globally unique
161  labelList procRegionOffset(Pstream::nProcs(), Zero);
162  for (label proci = 1; proci < Pstream::nProcs(); ++proci)
163  {
164  procRegionOffset[proci] =
165  procRegionOffset[proci - 1]
166  + procNRegions[proci - 1];
167  }
168 
169  // Apply the offset
170  for (label& regioni : newSetFaceRegions)
171  {
172  regioni += procRegionOffset[Pstream::myProcNo()];
173  }
174 
175  // Store the total number of regions across all processors
176  nRegions = procRegionOffset.last() + procNRegions.last();
177  }
178 
179  // Step 2: Create a region map which combines regions which are
180  // connected across coupled interfaces
181  labelList regionMap(identity(nRegions));
182  {
183  // Put region labels on connected boundary edges and synchronise to
184  // create a list of all regions connected to a given edge
185  labelListList meshEdgeRegions(mesh_.nEdges(), labelList());
186  forAll(newSetFaces, newSetFacei)
187  {
188  const label facei = newSetFaces[newSetFacei];
189  const label regioni = newSetFaceRegions[newSetFacei];
190  for (const label edgei : mesh_.faceEdges()[facei])
191  {
192  meshEdgeRegions[edgei] = labelList(one{}, regioni);
193  }
194  }
196  (
197  mesh_,
198  meshEdgeRegions,
199  ListOps::appendEqOp<label>(),
200  labelList()
201  );
202 
203  // Combine edge regions to create a list of what regions a given
204  // region is connected to
205  List<bitSet> regionRegions(nRegions);
206  forAll(newSetFaces, newSetFacei)
207  {
208  const label facei = newSetFaces[newSetFacei];
209  const label regioni = newSetFaceRegions[newSetFacei];
210  for (const label edgei : mesh_.faceEdges()[facei])
211  {
212  // Includes self region (removed below)
213  regionRegions[regioni].set(meshEdgeRegions[edgei]);
214  }
215  forAll(regionRegions, regioni)
216  {
217  // Remove self region
218  regionRegions[regioni].unset(regioni);
219  }
220  }
221  Pstream::listCombineGather(regionRegions, bitOrEqOp<bitSet>());
222  Pstream::listCombineScatter(regionRegions);
223 
224  // Collapse the region connections into a map between each region
225  // and the lowest numbered region that it connects to
226  forAll(regionRegions, regioni)
227  {
228  for (const label regi : regionRegions[regioni])
229  {
230  // minEqOp<label>()
231  regionMap[regi] = min(regionMap[regi], regionMap[regioni]);
232  }
233  }
234  }
235 
236  // Step 3: Combine connected regions
237  labelList regionNFaces;
238  {
239  // Remove duplicates from the region map
240  label regioni0 = 0;
241  forAll(regionMap, regioni)
242  {
243  if (regionMap[regioni] > regioni0)
244  {
245  ++ regioni0;
246  regionMap[regioni] = regioni0;
247  }
248  }
249 
250  // Recompute the number of regions
251  nRegions = regioni0 + 1;
252 
253  // Renumber the face region ID-s
254  newSetFaceRegions =
255  IndirectList<label>(regionMap, newSetFaceRegions);
256 
257  // Report the final number and size of the regions
258  regionNFaces = labelList(nRegions, Zero);
259  for (const label regioni : newSetFaceRegions)
260  {
261  ++ regionNFaces[regioni];
262  }
263  Pstream::listCombineGather(regionNFaces, plusEqOp<label>());
264  Pstream::listCombineScatter(regionNFaces);
265  Info<< " Found " << nRegions << " contiguous regions with "
266  << regionNFaces << " faces" << endl;
267  }
268 
269  // Step 4: Choose the closest region to output
270  label selectedRegioni = -1;
271  {
272  // Compute the region centres
273  scalarField regionWeights(nRegions, Zero);
274  pointField regionCentres(nRegions, Zero);
275  forAll(newSetFaces, newSetFacei)
276  {
277  const label facei = newSetFaces[newSetFacei];
278  const label regioni = newSetFaceRegions[newSetFacei];
279 
280  const scalar w = mag(mesh_.faceAreas()[facei]);
281  const point& c = mesh_.faceCentres()[facei];
282 
283  regionWeights[regioni] += w;
284  regionCentres[regioni] += w*c;
285  }
286  Pstream::listCombineGather(regionWeights, plusEqOp<scalar>());
287  Pstream::listCombineGather(regionCentres, plusEqOp<point>());
288  Pstream::listCombineScatter(regionWeights);
289  Pstream::listCombineScatter(regionCentres);
290  regionCentres /= regionWeights;
291 
292  // Find the region centroid closest to the reference point
293  selectedRegioni =
295  (
296  findMin(mag(regionCentres - point_)()),
297  minOp<label>()
298  );
299 
300  // Report the selection
301  Info<< " Selecting region " << selectedRegioni << " with "
302  << regionNFaces[selectedRegioni]
303  << " faces as the closest to point " << point_ << endl;
304  }
305 
306  // Step 5: Remove any faces from the set list not in the selected region
307  {
308  // Remove faces from the list by shuffling up and resizing
309  label newSetFacei0 = 0;
310  forAll(newSetFaces, newSetFacei)
311  {
312  newSetFaces[newSetFacei0] = newSetFaces[newSetFacei];
313 
314  if (newSetFaceRegions[newSetFacei] == selectedRegioni)
315  {
316  ++ newSetFacei0;
317  }
318  }
319  newSetFaces.resize(newSetFacei0);
320  }
321  }
322 
323  // Modify the face zone set
324  DynamicList<label> newAddressing;
325  DynamicList<bool> newFlipMap;
326  if (add)
327  {
328  // Start from copy
329  newAddressing = fzSet.addressing();
330  newFlipMap = fzSet.flipMap();
331 
332  // Add anything from the new set that is not already in the zone set
333  const auto& exclude = fzSet;
334  for (const label facei : newSetFaces)
335  {
336  if (!exclude.found(facei))
337  {
338  newAddressing.append(facei);
339  newFlipMap.append(cellIsAbovePlane[mesh_.faceOwner()[facei]]);
340  }
341  }
342  }
343  else
344  {
345  // Start from empty
346  newAddressing.reserve(fzSet.addressing().size());
347  newFlipMap.reserve(newAddressing.capacity());
348 
349  // Add everything from the zone set that is not also in the new set
350  bitSet exclude(newSetFaces);
351  for (const label facei : fzSet.addressing())
352  {
353  if (!exclude.found(facei))
354  {
355  newAddressing.append(facei);
356  newFlipMap.append(cellIsAbovePlane[mesh_.faceOwner()[facei]]);
357  }
358  }
359  }
360  fzSet.addressing().transfer(newAddressing);
361  fzSet.flipMap().transfer(newFlipMap);
362  fzSet.updateSet();
363 }
364 
365 
366 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
367 
369 (
370  const polyMesh& mesh,
371  const point& basePoint,
372  const vector& normal,
373  const faceAction action
374 )
375 :
377  point_(basePoint),
378  normal_(normal),
379  option_(action)
380 {}
381 
382 
384 (
385  const polyMesh& mesh,
386  const dictionary& dict
387 )
388 :
390  (
391  mesh,
392  dict.get<vector>("point"),
393  dict.get<vector>("normal"),
394  faceActionNames_.getOrDefault("option", dict, faceAction::ALL)
395  )
396 {}
397 
398 
400 (
401  const polyMesh& mesh,
402  Istream& is
403 )
404 :
406  point_(checkIs(is)),
407  normal_(checkIs(is)),
408  option_(faceActionNames_.read(checkIs(is)))
409 {}
410 
411 
412 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
413 
415 (
416  const topoSetSource::setAction action,
417  topoSet& set
418 ) const
419 {
420  if (!isA<faceZoneSet>(set))
421  {
423  << "Operation only allowed on a faceZoneSet." << endl;
424  return;
425  }
426 
427  faceZoneSet& zoneSet = refCast<faceZoneSet>(set);
428 
429  if (action == topoSetSource::NEW || action == topoSetSource::ADD)
430  {
431  if (verbose_)
432  {
433  Info<< " Adding faces that form a plane at "
434  << point_ << " with normal " << normal_ << endl;
435  }
436 
437  combine(zoneSet, true);
438  }
439  else if (action == topoSetSource::SUBTRACT)
440  {
441  if (verbose_)
442  {
443  Info<< " Removing faces that form a plane at "
444  << point_ << " with normal " << normal_ << endl;
445  }
446 
447  combine(zoneSet, false);
448  }
449 }
450 
451 
452 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::syncTools::syncEdgeList
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
Definition: syncToolsTemplates.C:800
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::planeToFaceZone::faceAction
faceAction
Enumeration defining the valid options.
Definition: planeToFaceZone.H:181
Foam::topoSetSource::ADD
Add elements to current set.
Definition: topoSetSource.H:103
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
uindirectPrimitivePatch.H
planeToFaceZone.H
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
PatchTools.H
faceZoneSet.H
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:215
Foam::topoSetSource::addToUsageTable
Class with constructor to add usage string to table.
Definition: topoSetSource.H:129
Foam::topoSetFaceZoneSource
The topoSetFaceZoneSource is a intermediate class for handling topoSet sources for selecting face zon...
Definition: topoSetFaceZoneSource.H:57
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::primitiveMesh::nEdges
label nEdges() const
Number of mesh edges.
Definition: primitiveMeshI.H:67
Foam::planeToFaceZone
A topoSetSource to select faces based on the adjacent cell centres spanning a given plane....
Definition: planeToFaceZone.H:174
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::topoSetSource::setAction
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:100
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::topoSetSource::NEW
Create a new set and ADD elements to it.
Definition: topoSetSource.H:104
polyMesh.H
Foam::planeToFaceZone::planeToFaceZone
planeToFaceZone()=delete
No default construct.
syncTools.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:396
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::ListListOps::combine
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
Definition: ListListOps.C:69
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::primitiveMesh::faceEdges
const labelListList & faceEdges() const
Definition: primitiveMeshEdges.C:528
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
Foam::faceZoneSet
Like faceSet but -reads data from faceZone -updates faceZone when writing.
Definition: faceZoneSet.H:51
Foam::topoSet
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:63
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::add
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Definition: FieldFieldFunctions.C:939
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::topoSetSource::SUBTRACT
Subtract elements from current set.
Definition: topoSetSource.H:105
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::planeToFaceZone::applyToSet
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
Apply specified action to the topoSet.
Definition: planeToFaceZone.C:415
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:77
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::PatchTools::markZones
static label markZones(const PrimitivePatch< FaceList, PointField > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
Foam::uindirectPrimitivePatch
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
Definition: uindirectPrimitivePatch.H:49
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::topoSetSource::mesh_
const polyMesh & mesh_
Reference to the mesh.
Definition: topoSetSource.H:156
Foam::findMin
label findMin(const ListType &input, label start=0)
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:89