targetVolumeToCell.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) 2012-2021 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 "targetVolumeToCell.H"
27 #include "polyMesh.H"
28 #include "globalMeshData.H"
29 #include "plane.H"
30 #include "cellSet.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::scalar Foam::targetVolumeToCell::volumeOfSet
45 (
46  const PackedBoolList& selected
47 ) const
48 {
49  scalar sumVol = 0.0;
50  forAll(selected, celli)
51  {
52  if (selected[celli])
53  {
54  sumVol += mesh_.cellVolumes()[celli];
55  }
56  }
57  return returnReduce(sumVol, sumOp<scalar>());
58 }
59 
60 
61 Foam::label Foam::targetVolumeToCell::selectCells
62 (
63  const scalar normalComp,
64  const PackedBoolList& maskSet,
65  PackedBoolList& selected
66 ) const
67 {
68  selected.setSize(mesh_.nCells());
69  selected = false;
70 
71  label nSelected = 0;
72 
73  forAll(mesh_.cellCentres(), celli)
74  {
75  const point& cc = mesh_.cellCentres()[celli];
76 
77  if (maskSet[celli] && ((cc&n_) < normalComp))
78  {
79  selected[celli] = true;
80  nSelected++;
81  }
82  }
83  return returnReduce(nSelected, sumOp<label>());
84 }
85 
86 
87 void Foam::targetVolumeToCell::combine(topoSet& set, const bool add) const
88 {
89  if (vol_ <= 0)
90  {
91  // Select no cells
92  return;
93  }
94 
95 
96  PackedBoolList maskSet(mesh_.nCells(), 1);
97  label nTotCells = mesh_.globalData().nTotalCells();
98  if (maskSetName_.size())
99  {
100  // Read cellSet
101  Info<< " Operating on subset defined by cellSet " << maskSetName_
102  << endl;
103 
104  maskSet = 0;
105  cellSet subset(mesh_, maskSetName_);
106 
107  forAllConstIter(cellSet, subset, iter)
108  {
109  maskSet[iter.key()] = 1;
110  }
111 
112  nTotCells = returnReduce(subset.size(), sumOp<label>());
113  }
114 
115 
116  // Get plane for min,max volume.
117  // Planes all have base (0 0 0) and fixed normal so work only on normal
118  // component.
119 
120  scalar maxComp = -great;
121  label maxCells = 0;
122  // scalar maxVol = 0;
123  scalar minComp = great;
124  {
125  const boundBox& bb = mesh_.bounds();
126  pointField points(bb.points());
127 
128  // label minPointi = -1;
129  label maxPointi = -1;
130  forAll(points, pointi)
131  {
132  scalar c = (points[pointi]&n_);
133  if (c > maxComp)
134  {
135  maxComp = c;
136  maxPointi = pointi;
137  }
138  else if (c < minComp)
139  {
140  minComp = c;
141  // minPointi = pointi;
142  }
143  }
144 
145  PackedBoolList maxSelected(mesh_.nCells());
146  maxCells = selectCells(maxComp, maskSet, maxSelected);
147  // maxVol = volumeOfSet(maxSelected);
148 
149  // Check that maxPoint indeed selects all cells
150  if (maxCells != nTotCells)
151  {
153  << "Plane " << plane(points[maxPointi], n_)
154  << " selects " << maxCells
155  << " cells instead of all " << nTotCells
156  << " cells. Results might be wrong." << endl;
157  }
158  }
159 
160 
161 
162  // Bisection
163  // ~~~~~~~~~
164 
165  PackedBoolList selected(mesh_.nCells());
166  label nSelected = -1;
167  scalar selectedVol = 0.0;
168  // scalar selectedComp = 0.0;
169 
170 
171  scalar low = minComp;
172  scalar high = maxComp;
173 
174  const scalar tolerance = small*100*(maxComp-minComp);
175 
176  while ((high-low) > tolerance)
177  {
178  scalar mid = 0.5*(low + high);
179 
180  nSelected = selectCells(mid, maskSet, selected);
181  selectedVol = volumeOfSet(selected);
182 
183  // Pout<< "High:" << high << " low:" << low << " mid:" << mid << nl
184  // << " nSelected:" << nSelected << nl
185  // << " vol :" << selectedVol << nl
186  // << endl;
187 
188  if (selectedVol < vol_)
189  {
190  low = mid;
191 
192  PackedBoolList highSelected(mesh_.nCells());
193  label nHigh = selectCells(high, maskSet, selected);
194  if (nSelected == nHigh)
195  {
196  break;
197  }
198  }
199  else
200  {
201  high = mid;
202 
203  PackedBoolList lowSelected(mesh_.nCells());
204  label nLow = selectCells(low, maskSet, selected);
205  if (nSelected == nLow)
206  {
207  break;
208  }
209  }
210  }
211 
212  nSelected = selectCells(high, maskSet, selected);
213  selectedVol = volumeOfSet(selected);
214 
215  if (selectedVol < vol_)
216  {
217  // selectedComp = high;
218  }
219  else
220  {
221  nSelected = selectCells(low, maskSet, selected);
222  selectedVol = volumeOfSet(selected);
223 
224  if (selectedVol < vol_)
225  {
226  // selectedComp = low;
227  }
228  else
229  {
231  << "Did not converge onto plane. " << nl
232  << "high plane:"
233  << plane(high*n_, n_)
234  << nl
235  << "low plane :"
236  << plane(low*n_, n_)
237  << endl;
238  }
239  }
240 
241 
242  Info<< " Selected " << nSelected << " with actual volume " << selectedVol
243  << endl;
244 
245  forAll(selected, celli)
246  {
247  if (selected[celli])
248  {
249  addOrDelete(set, celli, add);
250  }
251  }
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
256 
258 (
259  const polyMesh& mesh,
260  const scalar vol,
261  const vector& n
262 )
263 :
264  topoSetSource(mesh),
265  vol_(vol),
266  n_(n)
267 {}
268 
269 
271 (
272  const polyMesh& mesh,
273  const dictionary& dict
274 )
275 :
276  topoSetSource(mesh),
277  vol_(dict.lookup<scalar>("volume")),
278  n_(dict.lookup("normal")),
279  maskSetName_(dict.lookupOrDefault<word>("set", ""))
280 {}
281 
282 
283 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
284 
286 {}
287 
288 
289 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
290 
292 (
293  const topoSetSource::setAction action,
294  topoSet& set
295 ) const
296 {
297  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
298  {
299  Info<< " Adding cells up to target volume " << vol_
300  << " out of total volume " << gSum(mesh_.cellVolumes()) << endl;
301 
302  combine(set, true);
303  }
304  else if (action == topoSetSource::DELETE)
305  {
306  Info<< " Removing cells up to target volume " << vol_
307  << " out of total volume " << gSum(mesh_.cellVolumes()) << endl;
308 
309  combine(set, false);
310  }
311 }
312 
313 
314 // ************************************************************************* //
label n
#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.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const scalarField & cellVolumes() const
A topoSetSource to select cells based on the wanted volume of selected cells. Adapts a plane until it...
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
targetVolumeToCell(const polyMesh &mesh, const scalar vol, const vector &)
Construct from components.
virtual ~targetVolumeToCell()
Destructor.
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
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
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)
defineTypeNameAndDebug(combustionModel, 0)
treeBoundBox combine(const treeBoundBox &a, const treeBoundBox &b)
Definition: patchToPatch.C:78
static const char nl
Definition: Ostream.H:260
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dictionary dict