44 Foam::scalar Foam::targetVolumeToCell::volumeOfSet
46 const PackedBoolList& selected
63 const scalar normalComp,
64 const PackedBoolList& maskSet,
65 PackedBoolList& selected
68 selected.setSize(mesh_.nCells());
73 forAll(mesh_.cellCentres(), celli)
75 const point& cc = mesh_.cellCentres()[celli];
77 if (maskSet[celli] && ((cc&n_) < normalComp))
79 selected[celli] =
true;
87 void Foam::targetVolumeToCell::combine(topoSet& set,
const bool add)
const
96 PackedBoolList maskSet(mesh_.nCells(), 1);
97 label nTotCells = mesh_.globalData().nTotalCells();
98 if (maskSetName_.size())
101 Info<<
" Operating on subset defined by cellSet " << maskSetName_
105 cellSet
subset(mesh_, maskSetName_);
109 maskSet[iter.key()] = 1;
120 scalar maxComp = -great;
123 scalar minComp = great;
125 const boundBox& bb = mesh_.bounds();
129 label maxPointi = -1;
132 scalar
c = (
points[pointi]&n_);
138 else if (
c < minComp)
145 PackedBoolList maxSelected(mesh_.nCells());
146 maxCells = selectCells(maxComp, maskSet, maxSelected);
150 if (maxCells != nTotCells)
153 <<
"Plane " << plane(
points[maxPointi], n_)
154 <<
" selects " << maxCells
155 <<
" cells instead of all " << nTotCells
156 <<
" cells. Results might be wrong." <<
endl;
165 PackedBoolList selected(mesh_.nCells());
166 label nSelected = -1;
167 scalar selectedVol = 0.0;
171 scalar low = minComp;
172 scalar high = maxComp;
174 const scalar tolerance = small*100*(maxComp-minComp);
176 while ((high-low) > tolerance)
178 scalar mid = 0.5*(low + high);
180 nSelected = selectCells(mid, maskSet, selected);
181 selectedVol = volumeOfSet(selected);
188 if (selectedVol < vol_)
192 PackedBoolList highSelected(mesh_.nCells());
193 label nHigh = selectCells(high, maskSet, selected);
194 if (nSelected == nHigh)
203 PackedBoolList lowSelected(mesh_.nCells());
204 label nLow = selectCells(low, maskSet, selected);
205 if (nSelected == nLow)
212 nSelected = selectCells(high, maskSet, selected);
213 selectedVol = volumeOfSet(selected);
215 if (selectedVol < vol_)
221 nSelected = selectCells(low, maskSet, selected);
222 selectedVol = volumeOfSet(selected);
224 if (selectedVol < vol_)
231 <<
"Did not converge onto plane. " <<
nl
233 << plane(high*n_, n_)
242 Info<<
" Selected " << nSelected <<
" with actual volume " << selectedVol
249 addOrDelete(set, celli,
add);
277 vol_(
dict.lookup<scalar>(
"volume")),
278 n_(
dict.lookup(
"normal")),
279 maskSetName_(
dict.lookupOrDefault<
word>(
"set",
""))
299 Info<<
" Adding cells up to target volume " << vol_
300 <<
" out of total volume " <<
gSum(mesh_.cellVolumes()) <<
endl;
306 Info<<
" Removing cells up to target volume " << vol_
307 <<
" out of total volume " <<
gSum(mesh_.cellVolumes()) <<
endl;
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
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....
Mesh consisting of general polyhedral cells.
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.
setAction
Enumeration defining the valid actions.
General set of labels of mesh quantity (points, cells, faces).
A class for handling words, derived from string.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
vectorField pointField
pointField is a vectorField.
vector point
Point is a vector.
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)
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)