Go to the documentation of this file.
48 const treeBoundBox& bb,
52 OFstream str(time().
path()/
name +
".obj");
67 str <<
"l " <<
e[0] + 1 <<
' ' <<
e[1] + 1 <<
nl;
77 decomposition().overlapsThisProcessor(box)
78 || geometryToConformTo().overlaps(box);
81 return geometryToConformTo().overlaps(box);
90 decomposition().positionOnThisProcessor(
p)
91 && geometryToConformTo().inside(
p);
94 return geometryToConformTo().inside(
p);
106 return geometryToConformTo().wellInside
109 minimumSurfaceDistanceCoeffSqr_*
sqr(sizes)
113 Field<bool> inside(pts.size(),
true);
120 geometryToConformTo().wellInside
123 minimumSurfaceDistanceCoeffSqr_*
sqr(sizes)
129 decomposition().positionOnThisProcessor(pts)
148 inside[i] = (insideA[i] && insideB[i]);
165 inside = decomposition().positionOnThisProcessor(
p);
172 && geometryToConformTo().wellInside
175 minimumSurfaceDistanceCoeffSqr_*
sqr(size)
184 DynamicList<Vb::Point>& initialPoints,
185 const treeBoundBox& bb,
194 treeBoundBox subBB = bb.subBbox(i);
196 word newName = recursionName +
"_" +
Foam::name(i);
200 if (combinedOverlaps(subBB))
224 word(newName +
"_overlap")
227 Pout<< newName +
"_overlap " << subBB <<
endl;
230 if (!fillBox(initialPoints, subBB,
true))
247 else if (combinedInside(subBB.midpoint()))
257 Pout<< newName +
"_inside " << subBB <<
endl;
260 if (!fillBox(initialPoints, subBB,
false))
289 return treeDepth + 1;
295 DynamicList<Vb::Point>& initialPoints,
296 const treeBoundBox& bb,
300 const conformationSurfaces& geometry = geometryToConformTo();
302 label initialSize = initialPoints.size();
304 scalar maxCellSize = -GREAT;
306 scalar minCellSize = GREAT;
308 scalar maxDensity = 1/
pow3(minCellSize);
310 scalar volumeAdded = 0.0;
316 scalar totalVolume = bb.volume();
318 label trialPoints = 0;
320 bool wellInside =
false;
332 geometry.findSurfaceNearest
344 Pout<<
"box wellInside, no need to sample surface." <<
endl;
351 if (!overlapping && !wellInside)
360 scalarField cornerSizes = cellShapeControls().cellSize(corners);
362 Field<bool> insideCorners = combinedWellInside(corners, cornerSizes);
369 scalar
s = cornerSizes[i];
378 minCellSize =
max(
s, minCellSizeLimit_);
381 if (maxCellSize/minCellSize > maxSizeRatio_)
385 Pout<<
"Abort fill at corner sample stage,"
386 <<
" minCellSize " << minCellSize
387 <<
" maxCellSize " << maxCellSize
388 <<
" maxSizeRatio " << maxCellSize/minCellSize
395 if (!insideCorners[i])
402 Pout<<
"Inside box found to have some non-wellInside "
403 <<
"corners, using overlapping fill."
417 label nLine = 6*(surfRes_ - 2);
423 for (
label i = 0; i < surfRes_; i++)
427 for (
label j = 1; j < surfRes_ - 1 ; j++)
437 delta.x()*(surfRes_ - 1),
451 delta.y()*(surfRes_ - 1),
465 delta.z()*(surfRes_ - 1)
469 lineSizes = cellShapeControls().cellSize(linePoints);
471 Field<bool> insideLines = combinedWellInside
480 scalar
s = lineSizes[i];
489 minCellSize =
max(
s, minCellSizeLimit_);
492 if (maxCellSize/minCellSize > maxSizeRatio_)
496 Pout<<
"Abort fill at surface sample stage, "
497 <<
" minCellSize " << minCellSize
498 <<
" maxCellSize " << maxCellSize
499 <<
" maxSizeRatio " << maxCellSize/minCellSize
514 Pout<<
"Inside box found to have some non-"
515 <<
"wellInside surface points, using "
516 <<
"overlapping fill."
534 volRes_*volRes_*volRes_,
542 for (
label i = 0; i < volRes_; i++)
544 for (
label j = 0; j < volRes_; j++)
556 *(i + 0.5 + 0.1*(
rndGen().scalar01() - 0.5)),
558 *(j + 0.5 + 0.1*(
rndGen().scalar01() - 0.5)),
560 *(
k + 0.5 + 0.1*(
rndGen().scalar01() - 0.5))
571 scalarField sampleSizes = cellShapeControls().cellSize(samplePoints);
587 scalar
s = sampleSizes[i];
596 minCellSize =
max(
s, minCellSizeLimit_);
599 if (maxCellSize/minCellSize > maxSizeRatio_)
603 Pout<<
"Abort fill at sample stage,"
604 <<
" minCellSize " << minCellSize
605 <<
" maxCellSize " << maxCellSize
606 <<
" maxSizeRatio " << maxCellSize/minCellSize
619 Pout<<
"No sample points found inside box" <<
endl;
627 Pout<< scalar(nInside)/scalar(samplePoints.size())
628 <<
" full overlapping box" <<
endl;
631 totalVolume *= scalar(nInside)/scalar(samplePoints.size());
635 Pout<<
"Total volume to fill = " << totalVolume <<
endl;
642 maxDensity = 1/
pow3(
max(minCellSize, SMALL));
650 const point&
p = samplePoints[i];
652 scalar localSize = sampleSizes[i];
654 scalar localDensity = 1/
pow3(localSize);
665 if (localDensity/maxDensity >
rndGen().scalar01())
667 scalar localVolume = 1/localDensity;
669 if (volumeAdded + localVolume > totalVolume)
675 scalar addProbability =
676 (totalVolume - volumeAdded)/localVolume;
682 Pout<<
"totalVolume " << totalVolume <<
nl
683 <<
"volumeAdded " << volumeAdded <<
nl
684 <<
"localVolume " << localVolume <<
nl
685 <<
"addProbability " << addProbability <<
nl
690 if (addProbability > r)
703 volumeAdded += localVolume;
711 volumeAdded += localVolume;
717 if (volumeAdded < totalVolume)
721 Pout<<
"Adding random points, remaining volume "
722 << totalVolume - volumeAdded
726 maxDensity = 1/
pow3(
max(minCellSize, SMALL));
734 scalar localSize = cellShapeControls().cellSize(
p);
736 bool insidePoint =
false;
745 insidePoint = combinedWellInside(
p, localSize);
750 if (localSize > maxCellSize)
752 maxCellSize = localSize;
755 if (localSize < minCellSize)
757 minCellSize =
max(localSize, minCellSizeLimit_);
759 localSize = minCellSize;
763 maxDensity = 1/
pow3(
max(minCellSize, SMALL));
766 if (maxCellSize/minCellSize > maxSizeRatio_)
770 Pout<<
"Abort fill at random fill stage,"
771 <<
" minCellSize " << minCellSize
772 <<
" maxCellSize " << maxCellSize
773 <<
" maxSizeRatio " << maxCellSize/minCellSize
779 initialPoints.resize(initialSize);
784 scalar localDensity = 1/
pow3(
max(localSize, SMALL));
788 if (localDensity/maxDensity >
rndGen().scalar01())
790 scalar localVolume = 1/localDensity;
792 if (volumeAdded + localVolume > totalVolume)
798 scalar addProbability =
799 (totalVolume - volumeAdded)/localVolume;
805 Pout<<
"totalVolume " << totalVolume <<
nl
806 <<
"volumeAdded " << volumeAdded <<
nl
807 <<
"localVolume " << localVolume <<
nl
808 <<
"addProbability " << addProbability <<
nl
813 if (addProbability > r)
826 volumeAdded += localVolume;
834 volumeAdded += localVolume;
840 globalTrialPoints_ += trialPoints;
845 <<
" locations queried, " << initialPoints.size() - initialSize
846 <<
" points placed, ("
847 << scalar(initialPoints.size() - initialSize)
848 /scalar(
max(trialPoints, 1))
849 <<
" success rate)." <<
nl
850 <<
"minCellSize " << minCellSize
851 <<
", maxCellSize " << maxCellSize
852 <<
", ratio " << maxCellSize/minCellSize
864 const dictionary& initialPointsDict,
867 const conformationSurfaces& geometryToConformTo,
868 const cellShapeControl& cellShapeControls,
869 const autoPtr<backgroundMeshDecomposition>& decomposition
882 globalTrialPoints_(0),
885 detailsDict().lookupOrDefault<scalar>(
"minCellSizeLimit", 0.0)
892 detailsDict().lookupOrDefault<
label>(
"surfaceSampleResolution", volRes_)
895 if (maxSizeRatio_ <= 1.0)
900 <<
"The maxSizeRatio must be greater than one to be sensible, "
901 <<
"setting to " << maxSizeRatio_
935 Pout<<
" Filling box " << hierBB <<
endl;
952 reduce(nInitialPoints, sumOp<label>());
957 <<
indent << nInitialPoints <<
" points placed" <<
nl
961 <<
" success rate" <<
nl
964 <<
" levels of recursion (maximum)"
vectorField pointField
pointField is a vectorField.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const backgroundMeshDecomposition & decomposition() const
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
bool fillBox(DynamicList< Vb::Point > &initialPoints, const treeBoundBox &bb, bool overlapping) const
Fill the given box, optionally filling surface overlapping boxes.
#define forAll(list, i)
Loop across all elements in list.
void writeOBJ(Ostream &os, label &vertI, const tetPoints &tet)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static bool & parRun()
Is this a parallel run?
virtual List< Vb::Point > initialPoints() const
Return the initial points for the conformalVoronoiMesh.
void writeOBJ(const treeBoundBox &bb, fileName name) const
Write boundBox as obj.
Ostream & endl(Ostream &os)
Add newline and flush stream.
static const edgeList edges
Edge to point addressing.
void shuffle(UList< T > &)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
bool combinedOverlaps(const treeBoundBox &box) const
Check if the given box overlaps the geometry or, in parallel, the.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const conformationSurfaces & geometryToConformTo() const
Pre-declare SubField and related Field type.
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual const fileName & name() const
Return the name of the stream.
label minLevels_
Minimum normal level of recursion, can be more if a high density.
Macros for easy insertion into run-time selection tables.
scalar scalar01()
Returns the current sample.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
const double e
Elementary charge.
Ostream & indent(Ostream &os)
Indent stream.
Vector< scalar > vector
A scalar version of the templated Vector.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
label globalTrialPoints_
Trial points attempted to be placed in all boxes.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
prefixOSstream Pout(cout, "Pout")
const treeBoundBox & procBounds() const
Return the boundBox of this processor.
label k
Boltzmann constant.
bool combinedInside(const point &p) const
Check if the given point is inside the geometry and, in parallel,.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
label readLabel(Istream &is)
vector point
Point is a vector.
autoDensity(const dictionary &initialPointsDict, const Time &runTime, Random &rndGen, const conformationSurfaces &geometryToConformTo, const cellShapeControl &cellShapeControls, const autoPtr< backgroundMeshDecomposition > &decomposition)
Construct from components.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
cachedRandom rndGen(label(0), -1)
word name(const complex &)
Return a string representation of a complex.
Field< bool > combinedWellInside(const pointField &pts, const scalarField &sizes) const
Check if the given points are wellInside the geometry and, in.
stressControl lookup("compactNormalStress") >> compactNormalStress
label recurseAndFill(DynamicList< Vb::Point > &initialPoints, const treeBoundBox &bb, label levelLimit, word recursionName) const
Descend into octants of the supplied bound box.