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