targetVolumeToCell.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) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2018-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 "targetVolumeToCell.H"
30 #include "polyMesh.H"
31 #include "globalMeshData.H"
32 #include "plane.H"
33 #include "bitSet.H"
34 #include "cellSet.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(targetVolumeToCell, 0);
42  addToRunTimeSelectionTable(topoSetSource, targetVolumeToCell, word);
43  addToRunTimeSelectionTable(topoSetSource, targetVolumeToCell, istream);
44  addToRunTimeSelectionTable(topoSetCellSource, targetVolumeToCell, word);
45  addToRunTimeSelectionTable(topoSetCellSource, targetVolumeToCell, istream);
47  (
48  topoSetCellSource,
49  targetVolumeToCell,
50  word,
51  targetVolume
52  );
54  (
55  topoSetCellSource,
56  targetVolumeToCell,
57  istream,
58  targetVolume
59  );
60 }
61 
62 
63 Foam::topoSetSource::addToUsageTable Foam::targetVolumeToCell::usage_
64 (
65  targetVolumeToCell::typeName,
66  "\n Usage: targetVolumeToCell (nx ny nz)\n\n"
67  " Adjust plane until obtained selected volume\n\n"
68 );
69 
70 
71 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
72 
73 Foam::scalar Foam::targetVolumeToCell::volumeOfSet
74 (
75  const bitSet& selected
76 ) const
77 {
78  scalar sumVol = 0.0;
79 
80  // Loop over selected cells only
81  for (const label celli : selected)
82  {
83  sumVol += mesh_.cellVolumes()[celli];
84  }
85 
86  return returnReduce(sumVol, sumOp<scalar>());
87 }
88 
89 
90 Foam::label Foam::targetVolumeToCell::selectCells
91 (
92  const scalar normalComp,
93  const bitSet& maskSet,
94  bitSet& selected
95 ) const
96 {
97  selected.resize(mesh_.nCells());
98  selected = false;
99 
100  label nSelected = 0;
101 
102  forAll(mesh_.cellCentres(), celli)
103  {
104  const point& cc = mesh_.cellCentres()[celli];
105 
106  if (maskSet.test(celli) && ((cc & normal_) < normalComp))
107  {
108  selected.set(celli);
109  ++nSelected;
110  }
111  }
112 
113  return returnReduce(nSelected, sumOp<label>());
114 }
115 
116 
117 void Foam::targetVolumeToCell::combine(topoSet& set, const bool add) const
118 {
119  if (vol_ <= 0)
120  {
121  // Select no cells
122  return;
123  }
124 
125  bitSet maskSet(mesh_.nCells(), true);
126  label nTotCells = mesh_.globalData().nTotalCells();
127  if (maskSetName_.size())
128  {
129  // Read cellSet
130  if (verbose_)
131  {
132  Info<< " Operating on subset defined by cellSet "
133  << maskSetName_ << endl;
134  }
135 
136  maskSet = false;
137  cellSet subset(mesh_, maskSetName_);
138 
139  const labelHashSet& cellLabels = subset;
140  maskSet.setMany(cellLabels.begin(), cellLabels.end());
141 
142  nTotCells = returnReduce(subset.size(), sumOp<label>());
143  }
144 
145 
146  // Get plane for min,max volume.
147  // Planes all have base (0 0 0) and fixed normal so work only on normal
148  // component.
149 
150  scalar maxComp = -GREAT;
151  label maxCells = 0;
152  //scalar maxVol = 0;
153  scalar minComp = GREAT;
154  {
155  const boundBox& bb = mesh_.bounds();
156  pointField points(bb.points());
157 
158  //label minPointi = -1;
159  label maxPointi = -1;
160  forAll(points, pointi)
161  {
162  const scalar c = (points[pointi] & normal_);
163  if (c > maxComp)
164  {
165  maxComp = c;
166  maxPointi = pointi;
167  }
168  else if (c < minComp)
169  {
170  minComp = c;
171  //minPointi = pointi;
172  }
173  }
174 
175  bitSet maxSelected(mesh_.nCells());
176  maxCells = selectCells(maxComp, maskSet, maxSelected);
177  //maxVol = volumeOfSet(maxSelected);
178 
179  // Check that maxPoint indeed selects all cells
180  if (maxCells != nTotCells)
181  {
183  << "Plane " << plane(points[maxPointi], normal_)
184  << " selects " << maxCells
185  << " cells instead of all " << nTotCells
186  << " cells. Results might be wrong." << endl;
187  }
188  }
189 
190 
191  // Bisection
192  // ~~~~~~~~~
193 
194  bitSet selected(mesh_.nCells());
195  label nSelected = -1;
196  scalar selectedVol = 0.0;
197  //scalar selectedComp = 0.0;
198 
199 
200  scalar low = minComp;
201  scalar high = maxComp;
202 
203  const scalar tolerance = SMALL*100*(maxComp-minComp);
204 
205  while ((high-low) > tolerance)
206  {
207  const scalar mid = 0.5*(low + high);
208 
209  nSelected = selectCells(mid, maskSet, selected);
210  selectedVol = volumeOfSet(selected);
211 
212  //Pout<< "High:" << high << " low:" << low << " mid:" << mid << nl
213  // << " nSelected:" << nSelected << nl
214  // << " vol :" << selectedVol << nl
215  // << endl;
216 
217  if (selectedVol < vol_)
218  {
219  low = mid;
220 
221  bitSet highSelected(mesh_.nCells());
222  label nHigh = selectCells(high, maskSet, selected);
223  if (nSelected == nHigh)
224  {
225  break;
226  }
227  }
228  else
229  {
230  high = mid;
231 
232  bitSet lowSelected(mesh_.nCells());
233  label nLow = selectCells(low, maskSet, selected);
234  if (nSelected == nLow)
235  {
236  break;
237  }
238  }
239  }
240 
241  nSelected = selectCells(high, maskSet, selected);
242  selectedVol = volumeOfSet(selected);
243 
244  if (selectedVol < vol_)
245  {
246  //selectedComp = high;
247  }
248  else
249  {
250  nSelected = selectCells(low, maskSet, selected);
251  selectedVol = volumeOfSet(selected);
252 
253  if (selectedVol < vol_)
254  {
255  //selectedComp = low;
256  }
257  else
258  {
260  << "Did not converge onto plane. " << nl
261  << "high plane:"
262  << plane(high*normal_, normal_)
263  << nl
264  << "low plane :"
265  << plane(low*normal_, normal_)
266  << endl;
267  }
268  }
269 
270 
271  if (verbose_)
272  {
273  Info<< " Selected " << nSelected << " with actual volume "
274  << selectedVol << endl;
275  }
276 
277  addOrDelete(set, selected, add);
278 }
279 
280 
281 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
282 
284 (
285  const polyMesh& mesh,
286  const scalar vol,
287  const vector& normal,
288  const word& maskSetName
289 )
290 :
292  vol_(vol),
293  normal_(normal),
294  maskSetName_(maskSetName)
295 {}
296 
297 
299 (
300  const polyMesh& mesh,
301  const dictionary& dict
302 )
303 :
305  (
306  mesh,
307  dict.getCheck<scalar>("volume", scalarMinMax::ge(0)),
308  dict.get<vector>("normal"),
309  dict.getOrDefault<word>("set", "")
310  )
311 {}
312 
313 
315 (
316  const polyMesh& mesh,
317  Istream& is
318 )
319 :
321  vol_(readScalar(checkIs(is))),
322  normal_(checkIs(is))
323 {}
324 
325 
326 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
327 
329 (
330  const topoSetSource::setAction action,
331  topoSet& set
332 ) const
333 {
334  if (action == topoSetSource::ADD || action == topoSetSource::NEW)
335  {
336  if (verbose_)
337  {
338  Info<< " Adding cells up to target volume " << vol_
339  << " out of total volume "
340  << gSum(mesh_.cellVolumes()) << endl;
341  }
342 
343  combine(set, true);
344  }
345  else if (action == topoSetSource::SUBTRACT)
346  {
347  if (verbose_)
348  {
349  Info<< " Removing cells up to target volume " << vol_
350  << " out of total volume "
351  << gSum(mesh_.cellVolumes()) << endl;
352  }
353 
354  combine(set, false);
355  }
356 }
357 
358 
359 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::topoSetSource::ADD
Add elements to the set.
Definition: topoSetSource.H:101
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
globalMeshData.H
Foam::topoSetSource::addToUsageTable
Class with constructor to add usage string to table.
Definition: topoSetSource.H:124
Foam::subset
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
Foam::targetVolumeToCell::applyToSet
virtual void applyToSet(const topoSetSource::setAction action, topoSet &set) const
Apply specified action to the topoSet.
Definition: targetVolumeToCell.C:329
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::topoSetSource::setAction
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:99
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::topoSetSource::NEW
Create a new set and ADD elements to it.
Definition: topoSetSource.H:106
polyMesh.H
bitSet.H
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::targetVolumeToCell
A topoSetCellSource to select cells based on a target volume of cells. Adapts a plane until it has en...
Definition: targetVolumeToCell.H:159
plane.H
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 (uses stdout - output is on the master only)
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::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
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:121
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 the set.
Definition: topoSetSource.H:102
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::primitiveMesh::cellVolumes
const scalarField & cellVolumes() const
Definition: primitiveMeshCellCentresAndVols.C:96
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::MinMax::ge
static MinMax< T > ge(const T &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:31
Foam::Vector< scalar >
Foam::topoSetCellSource
The topoSetCellSource is a intermediate class for handling topoSet sources for selecting cells.
Definition: topoSetCellSource.H:54
points
const pointField & points
Definition: gmvOutputHeader.H:1
targetVolumeToCell.H
Foam::dictionary::getCheck
T getCheck(const word &keyword, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:94
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::targetVolumeToCell::targetVolumeToCell
targetVolumeToCell(const polyMesh &mesh, const scalar vol, const vector &normal, const word &maskSetName="")
Construct from components.
Definition: targetVolumeToCell.C:284
cellSet.H
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:409
Foam::topoSetSource::mesh_
const polyMesh & mesh_
Reference to the mesh.
Definition: topoSetSource.H:151
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303