Go to the documentation of this file.
50 targetVolumeToCell::typeName,
51 "\n Usage: targetVolumeToCell (nx ny nz)\n\n"
52 " Adjust plane until obtained selected volume\n\n"
68 sumVol += mesh_.cellVolumes()[cellI];
77 const scalar normalComp,
82 selected.
setSize(mesh_.nCells());
87 forAll(mesh_.cellCentres(), cellI)
89 const point& cc = mesh_.cellCentres()[cellI];
91 if (maskSet[cellI] && ((cc&n_) < normalComp))
93 selected[cellI] =
true;
123 maskSet[iter.key()] = 1;
134 scalar maxComp = -GREAT;
137 scalar minComp = GREAT;
143 label maxPointI = -1;
152 else if (
c < minComp)
160 maxCells =
selectCells(maxComp, maskSet, maxSelected);
164 if (maxCells != nTotCells)
168 <<
" selects " << maxCells
169 <<
" cells instead of all " << nTotCells
170 <<
" cells. Results might be wrong." <<
endl;
180 label nSelected = -1;
181 scalar selectedVol = 0.0;
185 scalar low = minComp;
186 scalar high = maxComp;
188 const scalar tolerance = SMALL*100*(maxComp-minComp);
190 while ((high-low) > tolerance)
192 scalar mid = 0.5*(low + high);
202 if (selectedVol <
vol_)
208 if (nSelected == nHigh)
219 if (nSelected == nLow)
229 if (selectedVol <
vol_)
238 if (selectedVol <
vol_)
245 <<
"Did not converge onto plane. " <<
nl
256 Info<<
" Selected " << nSelected <<
" with actual volume " << selectedVol
327 Info<<
" Adding cells up to target volume " << vol_
328 <<
" out of total volume " <<
gSum(mesh_.cellVolumes()) <<
endl;
334 Info<<
" Removing cells up to target volume " << vol_
335 <<
" out of total volume " <<
gSum(mesh_.cellVolumes()) <<
endl;
A class for handling words, derived from string.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
#define forAll(list, i)
Loop across all elements in list.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
Class with constructor to add usage string to table.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
label nTotalCells() const
Return total number of cells in decomposed mesh.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
label selectCells(const scalar normalComp, const PackedBoolList &, PackedBoolList &selected) const
Ostream & endl(Ostream &os)
Add newline and flush stream.
setAction
Enumeration defining the valid actions.
Type gSum(const FieldField< Field, Type > &f)
scalar volumeOfSet(const PackedBoolList &) const
const word maskSetName_
Optional name of cellSet to calculate volume in.
Mesh consisting of general polyhedral cells.
void setSize(const label, const unsigned int &val=0u)
Alias for resize()
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
const scalar vol_
Wanted volume.
virtual ~targetVolumeToCell()
Destructor.
General set of labels of mesh quantity (points, cells, faces).
A list of keyword definitions, which are a keyword followed by any number of values (e....
targetVolumeToCell(const polyMesh &mesh, const scalar vol, const vector &)
Construct from components.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Macros for easy insertion into run-time selection tables.
Base class of a source for a topoSet.
A collection of cell labels.
void combine(topoSet &set, const bool add) const
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
const boundBox & bounds() const
Return mesh bounding box.
const vector n_
Normal of plane to sweep.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
A bounding box defined in terms of the points at its extremities.
const dimensionedScalar c
Speed of light in a vacuum.
static addToUsageTable usage_
Add usage string.
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
const globalMeshData & globalData() const
Return parallel info.
tmp< pointField > points() const
Return corner points in an order corresponding to a 'hex' cell.
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
Report a warning using Foam::Warning.
void addOrDelete(topoSet &set, const label cellI, const bool) const
Add (if bool) cellI to set or delete cellI from set.