49void Foam::patchSeedSet::calcPatchSamples
51 const label nAvailable,
52 const label nPatchPoints,
53 DynamicList<point>& samplingPts,
54 DynamicList<label>& samplingCells,
55 DynamicList<label>& samplingFaces,
56 DynamicList<label>& samplingSegments,
57 DynamicList<scalar>& samplingCurveDist
65 Random&
rndGen = *rndGenPtr_;
67 globalIndex globalSampleNumbers(nAvailable);
68 label nGlobalPatchPoints =
returnReduce(nPatchPoints, sumOp<label>());
74 const bool perturb =
true;
76 for (
const label patchi : patchSet_)
79 triangulatedPatch tp(
pp, perturb);
81 const label np = nAvailable*
pp.size()/scalar(nGlobalPatchPoints);
82 for (label i = 0; i < np; ++i)
84 tp.randomLocalPoint(
rndGen, pt, facei, celli);
86 samplingPts.append(pt);
87 samplingCells.append(celli);
88 samplingFaces.append(
pp.start() + facei);
89 samplingSegments.append(0);
90 samplingCurveDist.append(globalSampleNumbers.toGlobal(i));
96void Foam::patchSeedSet::calcSelectedLocations
98 const label nAvailable,
99 const label nPatchPoints,
116 for (
const label patchi : patchSet_)
121 patchFaces[sz++] =
pp.start() + localFacei;
160 const scalar globalDistSqr
166 for (label sampleI = 0; sampleI < nAvailable; ++sampleI)
168 const auto& treeData = boundaryTree.shapes();
169 const point& sample = selectedLocations_[sampleI];
172 auto& distSqrProc = nearest[sampleI].second();
174 nearInfo = boundaryTree.findNearest
187 nearInfo.setPoint(treeData.centre(nearInfo.index()));
189 distSqrProc.first() = sample.distSqr(nearInfo.point());
204 if (nearest[sampleI].first().hit())
206 label procI = nearest[sampleI].second().second();
207 label index = nearest[sampleI].first().index();
219 Pout<<
"Found " << newPatchFaces.size()
220 <<
" out of " << nAvailable
221 <<
" on local processor" <<
endl;
224 patchFaces.transfer(newPatchFaces);
231 if (totalSize > nAvailable)
234 label myMaxPoints = scalar(patchFaces.size())/totalSize*nAvailable;
237 for (label iter = 0; iter < 4; ++iter)
247 subset.setSize(myMaxPoints);
253 Pout<<
"In random mode : selected " <<
subset.size()
254 <<
" faces out of " << patchFaces.size() <<
endl;
262 globalIndex globalSampleNumbers(patchFaces.size());
264 samplingPts.setCapacity(patchFaces.size());
265 samplingCells.setCapacity(patchFaces.size());
266 samplingFaces.setCapacity(patchFaces.size());
267 samplingSegments.setCapacity(patchFaces.size());
268 samplingCurveDist.setCapacity(patchFaces.size());
275 const label facei = patchFaces[i];
293 info.point() + 1
e-1*(cc-info.point())
298 samplingPts.append(info.point());
300 samplingCells.
append(celli);
301 samplingFaces.append(facei);
302 samplingSegments.append(0);
303 samplingCurveDist.append(globalSampleNumbers.toGlobal(i));
308void Foam::patchSeedSet::genSamples()
326 samplingPts.shrink();
327 samplingCells.shrink();
328 samplingFaces.shrink();
329 samplingSegments.shrink();
330 samplingCurveDist.shrink();
335 std::move(samplingPts),
336 std::move(samplingCells),
337 std::move(samplingFaces),
338 std::move(samplingSegments),
339 std::move(samplingCurveDist)
349void Foam::patchSeedSet::calcSamples
362 rndGenPtr_.reset(
new Random(0));
365 label nPatchPoints = 0;
366 for (
const label patchi : patchSet_)
369 nPatchPoints +=
pp.
size();
374 label nAvailable =
min(maxPoints_, selectedLocations_.size());
376 calcSelectedLocations
387 nAvailable = maxPoints_ - nAvailable;
417 maxPoints_(
dict.get<label>(
"maxPoints")),
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
const Addr & addressing() const noexcept
The addressing used for the list.
label size() const noexcept
The number of elements in the list.
A List with indirect addressing.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void append(const T &val)
Append an element at the end of the list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
static void listCombineReduce(UList< T > &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::listReduce with an in-place cop.
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
bool get(const label i) const
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const word & name() const noexcept
The coord-set name.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Non-pointer based hierarchical recursive searching.
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Initialises points on or just off patch.
patchSeedSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
virtual const labelList & faceOwner() const
Return face owner.
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
A patch is a list of labels that address the faces in the global face list.
const vectorField & cellCentres() const
Holds list of sampling points which is filled at construction time. Various implementations of this b...
const meshSearch & searchEngine() const noexcept
sampledSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const coordSet::coordFormat axisType)
Construct from components.
const polyMesh & mesh() const noexcept
Standard boundBox with extra functionality for use in octree.
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Encapsulation of data for searching on faces.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define DebugInfo
Report an information message using Foam::Info.
Namespace for handling debugging switches.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
UIndirectList< label > labelUIndList
UIndirectList of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
vector point
Point is a vector.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
#define forAll(list, i)
Loop across all elements in list.