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-2019 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 {
37  defineTypeNameAndDebug(targetVolumeToCell, 0);
38  addToRunTimeSelectionTable(topoSetSource, targetVolumeToCell, word);
39  addToRunTimeSelectionTable(topoSetSource, targetVolumeToCell, istream);
40 }
41 
42 
43 Foam::topoSetSource::addToUsageTable Foam::targetVolumeToCell::usage_
44 (
45  targetVolumeToCell::typeName,
46  "\n Usage: targetVolumeToCell (nx ny nz)\n\n"
47  " Adjust plane until obtained selected volume\n\n"
48 );
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 Foam::scalar Foam::targetVolumeToCell::volumeOfSet
54 (
55  const PackedBoolList& selected
56 ) const
57 {
58  scalar sumVol = 0.0;
59  forAll(selected, celli)
60  {
61  if (selected[celli])
62  {
63  sumVol += mesh_.cellVolumes()[celli];
64  }
65  }
66  return returnReduce(sumVol, sumOp<scalar>());
67 }
68 
69 
70 Foam::label Foam::targetVolumeToCell::selectCells
71 (
72  const scalar normalComp,
73  const PackedBoolList& maskSet,
74  PackedBoolList& selected
75 ) const
76 {
77  selected.setSize(mesh_.nCells());
78  selected = false;
79 
80  label nSelected = 0;
81 
82  forAll(mesh_.cellCentres(), celli)
83  {
84  const point& cc = mesh_.cellCentres()[celli];
85 
86  if (maskSet[celli] && ((cc&n_) < normalComp))
87  {
88  selected[celli] = true;
89  nSelected++;
90  }
91  }
92  return returnReduce(nSelected, sumOp<label>());
93 }
94 
95 
96 void Foam::targetVolumeToCell::combine(topoSet& set, const bool add) const
97 {
98  if (vol_ <= 0)
99  {
100  // Select no cells
101  return;
102  }
103 
104 
105  PackedBoolList maskSet(mesh_.nCells(), 1);
106  label nTotCells = mesh_.globalData().nTotalCells();
107  if (maskSetName_.size())
108  {
109  // Read cellSet
110  Info<< " Operating on subset defined by cellSet " << maskSetName_
111  << endl;
112 
113  maskSet = 0;
114  cellSet subset(mesh_, maskSetName_);
115 
116  forAllConstIter(cellSet, subset, iter)
117  {
118  maskSet[iter.key()] = 1;
119  }
120 
121  nTotCells = returnReduce(subset.size(), sumOp<label>());
122  }
123 
124 
125  // Get plane for min,max volume.
126  // Planes all have base (0 0 0) and fixed normal so work only on normal
127  // component.
128 
129  scalar maxComp = -great;
130  label maxCells = 0;
131  // scalar maxVol = 0;
132  scalar minComp = great;
133  {
134  const boundBox& bb = mesh_.bounds();
135  pointField points(bb.points());
136 
137  // label minPointi = -1;
138  label maxPointi = -1;
139  forAll(points, pointi)
140  {
141  scalar c = (points[pointi]&n_);
142  if (c > maxComp)
143  {
144  maxComp = c;
145  maxPointi = pointi;
146  }
147  else if (c < minComp)
148  {
149  minComp = c;
150  // minPointi = pointi;
151  }
152  }
153 
154  PackedBoolList maxSelected(mesh_.nCells());
155  maxCells = selectCells(maxComp, maskSet, maxSelected);
156  // maxVol = volumeOfSet(maxSelected);
157 
158  // Check that maxPoint indeed selects all cells
159  if (maxCells != nTotCells)
160  {
162  << "Plane " << plane(points[maxPointi], n_)
163  << " selects " << maxCells
164  << " cells instead of all " << nTotCells
165  << " cells. Results might be wrong." << endl;
166  }
167  }
168 
169 
170 
171  // Bisection
172  // ~~~~~~~~~
173 
174  PackedBoolList selected(mesh_.nCells());
175  label nSelected = -1;
176  scalar selectedVol = 0.0;
177  // scalar selectedComp = 0.0;
178 
179 
180  scalar low = minComp;
181  scalar high = maxComp;
182 
183  const scalar tolerance = small*100*(maxComp-minComp);
184 
185  while ((high-low) > tolerance)
186  {
187  scalar mid = 0.5*(low + high);
188 
189  nSelected = selectCells(mid, maskSet, selected);
190  selectedVol = volumeOfSet(selected);
191 
192  // Pout<< "High:" << high << " low:" << low << " mid:" << mid << nl
193  // << " nSelected:" << nSelected << nl
194  // << " vol :" << selectedVol << nl
195  // << endl;
196 
197  if (selectedVol < vol_)
198  {
199  low = mid;
200 
201  PackedBoolList highSelected(mesh_.nCells());
202  label nHigh = selectCells(high, maskSet, selected);
203  if (nSelected == nHigh)
204  {
205  break;
206  }
207  }
208  else
209  {
210  high = mid;
211 
212  PackedBoolList lowSelected(mesh_.nCells());
213  label nLow = selectCells(low, maskSet, selected);
214  if (nSelected == nLow)
215  {
216  break;
217  }
218  }
219  }
220 
221  nSelected = selectCells(high, maskSet, selected);
222  selectedVol = volumeOfSet(selected);
223 
224  if (selectedVol < vol_)
225  {
226  // selectedComp = high;
227  }
228  else
229  {
230  nSelected = selectCells(low, maskSet, selected);
231  selectedVol = volumeOfSet(selected);
232 
233  if (selectedVol < vol_)
234  {
235  // selectedComp = low;
236  }
237  else
238  {
240  << "Did not converge onto plane. " << nl
241  << "high plane:"
242  << plane(high*n_, n_)
243  << nl
244  << "low plane :"
245  << plane(low*n_, n_)
246  << endl;
247  }
248  }
249 
250 
251  Info<< " Selected " << nSelected << " with actual volume " << selectedVol
252  << endl;
253 
254  forAll(selected, celli)
255  {
256  if (selected[celli])
257  {
258  addOrDelete(set, celli, add);
259  }
260  }
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
265 
267 (
268  const polyMesh& mesh,
269  const scalar vol,
270  const vector& n
271 )
272 :
273  topoSetSource(mesh),
274  vol_(vol),
275  n_(n)
276 {}
277 
278 
280 (
281  const polyMesh& mesh,
282  const dictionary& dict
283 )
284 :
285  topoSetSource(mesh),
286  vol_(dict.lookup<scalar>("volume")),
287  n_(dict.lookup("normal")),
288  maskSetName_(dict.lookupOrDefault<word>("set", ""))
289 {}
290 
291 
293 (
294  const polyMesh& mesh,
295  Istream& is
296 )
297 :
298  topoSetSource(mesh),
299  vol_(readScalar(checkIs(is))),
300  n_(checkIs(is))
301 {}
302 
303 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
304 
306 {}
307 
308 
309 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
310 
312 (
313  const topoSetSource::setAction action,
314  topoSet& set
315 ) const
316 {
317  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
318  {
319  Info<< " Adding cells up to target volume " << vol_
320  << " out of total volume " << gSum(mesh_.cellVolumes()) << endl;
321 
322  combine(set, true);
323  }
324  else if (action == topoSetSource::DELETE)
325  {
326  Info<< " Removing cells up to target volume " << vol_
327  << " out of total volume " << gSum(mesh_.cellVolumes()) << endl;
328 
329  combine(set, false);
330  }
331 }
332 
333 
334 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual ~targetVolumeToCell()
Destructor.
const dimensionedScalar & c
Speed of light in a vacuum.
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
Definition: ListListOps.C:34
Macros for easy insertion into run-time selection tables.
Base class of a source for a topoSet.
Definition: topoSetSource.H:63
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Type gSum(const FieldField< Field, Type > &f)
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:82
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
Class with constructor to add usage string to table.
targetVolumeToCell(const polyMesh &mesh, const scalar vol, const vector &)
Construct from components.
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812