57void Foam::searchableSurfaceCollection::findNearest
60 scalarField& minDistSqr,
61 List<pointIndexHit>& nearestInfo,
62 labelList& nearestSurf
66 nearestInfo.setSize(
samples.size());
68 nearestSurf.setSize(
samples.size());
71 List<pointIndexHit> hitInfo(
samples.size());
77 subGeom_[surfI].findNearest
81 transform_[surfI].localPosition(
samples),
90 if (hitInfo[pointi].hit())
94 point globalPt = transform_[surfI].globalPosition
98 hitInfo[pointi].
point(),
103 scalar distSqr = globalPt.distSqr(
samples[pointi]);
105 if (distSqr < minDistSqr[pointi])
107 minDistSqr[pointi] = distSqr;
108 nearestInfo[pointi].hitPoint(globalPt);
109 nearestInfo[pointi].setIndex
111 hitInfo[pointi].index()
112 + indexOffset_[surfI]
114 nearestSurf[pointi] = surfI;
124void Foam::searchableSurfaceCollection::sortHits
136 if (info[pointi].hit())
138 label index = info[pointi].index();
139 label surfI =
findLower(indexOffset_, index+1);
145 surfInfo.setSize(subGeom_.size());
147 infoMap.setSize(subGeom_.size());
151 surfInfo[surfI].setSize(nHits[surfI]);
152 infoMap[surfI].setSize(nHits[surfI]);
158 if (info[pointi].hit())
160 label index = info[pointi].index();
161 label surfI =
findLower(indexOffset_, index+1);
165 label localI = nHits[surfI]++;
169 index-indexOffset_[surfI]
171 infoMap[surfI][localI] = pointi;
179Foam::searchableSurfaceCollection::searchableSurfaceCollection
186 instance_(
dict.size()),
188 transform_(
dict.size()),
189 subGeom_(
dict.size()),
190 mergeSubRegions_(
dict.
get<bool>(
"mergeSubRegions")),
191 indexOffset_(
dict.size()+1)
196 label startIndex = 0;
197 for (
const entry& dEntry :
dict)
201 instance_[surfI] = dEntry.keyword();
205 sDict.readEntry(
"scale", scale_[surfI]);
211 new coordSystem::cartesian(sDict,
"transform")
214 const word subGeomName(sDict.get<word>(
"surface"));
217 searchableSurface&
s =
218 io.db().lookupObjectRef<searchableSurface>(subGeomName);
223 if (
s.size() !=
s.globalSize())
226 <<
"Cannot use a distributed surface in a collection."
230 subGeom_.set(surfI, &
s);
232 indexOffset_[surfI] = startIndex;
233 startIndex += subGeom_[surfI].size();
235 Info<<
" instance : " << instance_[surfI] <<
endl;
236 Info<<
" surface : " <<
s.name() <<
endl;
237 Info<<
" scale : " << scale_[surfI] <<
endl;
238 Info<<
" transform: " << transform_[surfI] <<
endl;
243 indexOffset_[surfI] = startIndex;
245 instance_.setSize(surfI);
246 scale_.setSize(surfI);
247 transform_.setSize(surfI);
248 subGeom_.setSize(surfI);
249 indexOffset_.setSize(surfI+1);
256 const boundBox& surfBb = subGeom_[surfI].bounds();
259 const point surfBbMin = transform_[surfI].globalPosition
267 const point surfBbMax = transform_[surfI].globalPosition
292 if (regions_.empty())
294 regionOffset_.setSize(subGeom_.size());
299 regionOffset_[surfI] = allRegions.
size();
301 if (mergeSubRegions_)
308 const wordList& subRegions = subGeom_[surfI].regions();
324 return indexOffset_.last();
332 auto& ctrs = tctrs.ref();
343 ctrs[coordI++] = transform_[surfI].globalPosition
364 centres.setSize(size());
365 radiusSqr.setSize(centres.size());
372 scalar maxScale =
cmptMax(scale_[surfI]);
376 subGeom_[surfI].boundingSpheres(subCentres, subRadiusSqr);
380 centres[coordI] = transform_[surfI].globalPosition
388 radiusSqr[coordI] = maxScale*subRadiusSqr[i];
395Foam::tmp<Foam::pointField>
403 nPoints += subGeom_[surfI].points()().size();
407 auto&
pts = tpts.ref();
418 pts[
nPoints++] = transform_[surfI].globalPosition
433void Foam::searchableSurfaceCollection::findNearest
474 transform_[surfI].localPosition
484 transform_[surfI].localPosition
491 subGeom_[surfI].findLine(e0, e1, hitInfo);
495 if (hitInfo[pointi].hit())
498 nearest[pointi] = transform_[surfI].globalPosition
502 hitInfo[pointi].
point(),
506 info[pointi] = hitInfo[pointi];
507 info[pointi].point() = nearest[pointi];
508 info[pointi].setIndex
510 hitInfo[pointi].index()
511 + indexOffset_[surfI]
523 if (info[pointi].hit())
525 vector n(end[pointi] - start[pointi]);
526 scalar magN =
mag(
n);
532 scalar
s = ((info[pointi].point()-start[pointi])&
n);
537 <<
"point:" << info[pointi]
539 <<
" outside vector "
540 <<
" start:" << start[pointi]
541 <<
" end:" <<
end[pointi]
572 findLine(start, end, nearestInfo);
577 if (nearestInfo[pointi].hit())
580 info[pointi][0] = nearestInfo[pointi];
584 info[pointi].clear();
596 if (subGeom_.size() == 0)
598 else if (subGeom_.size() == 1)
600 if (mergeSubRegions_)
602 region.setSize(info.size());
603 region = regionOffset_[0];
607 subGeom_[0].getRegion(info, region);
618 sortHits(info, surfInfo, infoMap);
620 region.setSize(info.size());
625 if (mergeSubRegions_)
633 region[map[i]] = regionOffset_[surfI];
642 subGeom_[surfI].getRegion(surfInfo[surfI], surfRegion);
647 region[map[i]] = regionOffset_[surfI] + surfRegion[i];
661 if (subGeom_.size() == 0)
663 else if (subGeom_.size() == 1)
665 subGeom_[0].getNormal(info, normal);
675 sortHits(info, surfInfo, infoMap);
677 normal.setSize(info.size());
683 subGeom_[surfI].getNormal(surfInfo[surfI], surfNormal);
686 surfNormal = transform_[surfI].globalVector(surfNormal);
691 normal[map[i]] = surfNormal[i];
705 <<
"Volume type not supported for collection."
713 const bool keepNonLocal,
730 subGeom_[surfI].distribute
745 subGeom_[surfI].setField
752 subGeom_[surfI].size(),
763 const List<pointIndexHit>& info,
767 if (subGeom_.size() == 0)
769 else if (subGeom_.size() == 1)
771 subGeom_[0].getField(info, values);
781 sortHits(info, surfInfo, infoMap);
787 subGeom_[surfI].getField(surfInfo[surfI], surfValues);
789 if (surfValues.size())
792 values.setSize(info.size());
797 values[map[i]] = surfValues[i];
Various functions to operate on Lists.
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
void append(const T &val)
Copy append an element to the end of this list.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
void setSize(label n)
Alias for resize().
A non-owning sub-view of a List (allocated or unallocated storage).
void size(const label n)
Older name for setAddressableSize.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A bounding box defined in terms of min/max extrema points.
const point & max() const noexcept
Maximum describing the bounding box.
const point & min() const noexcept
Minimum describing the bounding box.
void reset()
Reset to an inverted box.
A Cartesian coordinate system.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
A keyword and a list of tokens is an 'entry'.
Set of transformed searchableSurfaces. Does not do boolean operations so when meshing might find part...
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
virtual void distribute(const List< treeBoundBox > &, const bool keepNonLocal, autoPtr< mapDistribute > &faceMap, autoPtr< mapDistribute > &pointMap)
Set bounds of surface. Bounds currently set as list of.
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &) const
Get all intersections in order from start to end.
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
virtual label size() const
Range of local indices that can be returned.
virtual void boundingSpheres(pointField ¢res, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
virtual void setField(const labelList &values)
WIP. Store element-wise field.
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
virtual void getRegion(const List< pointIndexHit > &, labelList ®ion) const
From a set of points and indices get the region.
virtual const wordList & regions() const
Names of regions.
virtual ~searchableSurfaceCollection()
Destructor.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
virtual tmp< pointField > points() const
Get the points that define the surface.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
virtual const boundBox & bounds() const
Return const reference to boundBox.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< word > wordList
List of word.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
messageStream Info
Information stream (stdout output on master, null elsewhere).
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define forAll(list, i)
Loop across all elements in list.
scalarField samples(nIntervals, Zero)