44Foam::sampledMeshedSurface::samplingSourceNames_
46 { samplingSource::cells,
"cells" },
47 { samplingSource::insideCells,
"insideCells" },
48 { samplingSource::boundaryFaces,
"boundaryFaces" },
100void Foam::sampledMeshedSurface::setZoneMap()
106 const auto& zones =
s.surfZones();
111 if (zoneIds_.
empty() || zones.size() <= 1)
122 const label len =
min(zones[zonei].
size(), zoneIds_.
size() - beg);
125 SubList<label>(zoneIds_, len, beg) = zonei;
132 const label len = (zoneIds_.
size() - beg);
136 SubList<label>(zoneIds_, len, beg) =
max(0, zones.size()-1);
148 globalIndex globalCells(onBoundary() ?
mesh().
nFaces() :
mesh().nCells());
152 const pointField& fc = surface_.faceCentres();
155 typedef Tuple2<scalar, label> nearInfo;
158 if (sampleSource_ == samplingSource::cells)
162 const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
166 const point& pt = fc[facei];
167 auto& near = nearest[facei];
173 near.first() =
magSqr(info.hitPoint()-pt);
174 near.second() = globalCells.toGlobal(info.index());
178 else if (sampleSource_ == samplingSource::insideCells)
182 const auto& cellTree = meshSearcher.cellTree();
186 const point& pt = fc[facei];
187 auto& near = nearest[facei];
189 if (cellTree.bb().contains(pt))
191 const label index = cellTree.findInside(pt);
195 near.second() = globalCells.toGlobal(index);
205 const auto& bndTree = meshSearcher.nonCoupledBoundaryTree();
209 const point& pt = fc[facei];
210 auto& near = nearest[facei];
216 near.first() =
magSqr(info.hitPoint()-pt);
220 bndTree.shapes().faceLabels()[info.index()]
232 labelList cellOrFaceLabels(fc.size(), -1);
234 bitSet facesToSubset(fc.size());
238 const auto& near = nearest[facei];
240 const label index = near.second();
247 else if (globalCells.isLocal(index))
249 facesToSubset.set(facei);
252 cellOrFaceLabels[facei] =
254 (near.first() < maxDistanceSqr_)
255 ? globalCells.toLocal(index)
264 Pout<<
"Local out of faces:" << cellOrFaceLabels.size()
265 <<
" keeping:" << facesToSubset.count() <<
endl;
275 s = surface_.subsetMesh(facesToSubset, pointMap,
faceMap);
294 sampleElements_.resize(pointMap.size(), -1);
297 labelList pointToFace(std::move(pointMap));
301 const face&
f =
s[facei];
303 for (
const label labi :
f)
305 pointToFace[labi] = facei;
310 if (sampleSource_ == samplingSource::cells)
315 forAll(samplePoints_, pointi)
318 const point pt = samplePoints_[pointi];
320 const label celli = cellOrFaceLabels[pointToFace[pointi]];
322 sampleElements_[pointi] = celli;
327 && !
mesh().pointInCell(pt, celli, meshSearcher.decompMode())
333 scalar minDistSqr = VGREAT;
335 for (
const label facei :
mesh().
cells()[celli])
346 if (info.distance() < minDistSqr)
348 minDistSqr = info.distance();
349 samplePoints_[pointi] = info.rawPoint();
355 else if (sampleSource_ == samplingSource::insideCells)
360 forAll(samplePoints_, pointi)
362 const label celli = cellOrFaceLabels[pointToFace[pointi]];
364 sampleElements_[pointi] = celli;
373 forAll(samplePoints_, pointi)
375 const point& pt = samplePoints_[pointi];
377 const label facei = cellOrFaceLabels[pointToFace[pointi]];
379 sampleElements_[pointi] = facei;
383 samplePoints_[pointi] =
405 sampleElements_.
transfer(cellOrFaceLabels);
406 samplePoints_.clear();
413 OFstream str(
mesh().time().
path()/
"surfaceToMesh.obj");
414 Info<<
"Dumping correspondence from local surface (points or faces)"
415 <<
" to mesh (cells or faces) to " << str.name() <<
endl;
427 forAll(samplePoints_, pointi)
435 label elemi = sampleElements_[pointi];
438 str <<
"l " << vertI-2 <<
' ' << vertI-1 <<
' ' << vertI <<
nl;
444 forAll(sampleElements_, triI)
449 label elemi = sampleElements_[triI];
452 str <<
"l " << vertI-1 <<
' ' << vertI <<
nl;
457 needsUpdate_ =
false;
468 const word& surfaceName,
474 surfaceName_(surfaceName),
480 sampleSource_(sampleSource),
486 maxDistanceSqr_(
Foam::
sqr(GREAT)),
513 sampleSource_(samplingSourceNames_.get(
"source",
dict)),
515 keepIds_(
dict.getOrDefault(
"keepIds", true)),
519 maxDistanceSqr_(
Foam::
sqr(GREAT)),
520 defaultValues_(
dict.subOrEmptyDict(
"defaultValue"))
534 maxDistanceSqr_ =
Foam::sqr(maxDistanceSqr_);
539 includePatches.
uniq();
544 if (!includePatches.
empty())
546 Info<<
"Subsetting surface " << surfaceName_
565 for (
const label zonei : zoneIndices)
571 if (includeMap.
none())
574 <<
"Patch selection results in an empty surface"
575 <<
" - ignoring" <<
nl;
577 else if (!includeMap.
all())
584 <<
"Bad surface subset (empty)"
585 <<
" - skip and hope for the best" <<
nl;
617 sampleElements_.clear();
618 samplePoints_.clear();
633 treeBoundBox bb(surface_.points(), surface_.meshPoints());
644 <<
"Surface " << surfaceName_
645 <<
" does not overlap bounding box of mesh " <<
mesh().
bounds()
651 bb.
min() += (0.5-1
e-6)*span;
652 bb.
max() -= (0.5-1
e-6)*span;
658 bb.
min() -= 0.5*span;
659 bb.
max() += 0.5*span;
667 return update(meshSearcher);
681 return update(meshSearcher);
690 return sampleOnFaces(sampler);
699 return sampleOnFaces(sampler);
708 return sampleOnFaces(sampler);
717 return sampleOnFaces(sampler);
726 return sampleOnFaces(sampler);
735 return sampleOnPoints(interpolator);
744 return sampleOnPoints(interpolator);
752 return sampleOnPoints(interpolator);
761 return sampleOnPoints(interpolator);
770 return sampleOnPoints(interpolator);
776 os <<
"meshedSurface: " <<
name() <<
" :"
777 <<
" surface:" << surfaceName_;
781 os <<
" faces:" << faces().size()
783 <<
" zoneids:" << zoneIds().size();
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Minimal example by using system/controlDict.functions:
Defines the attributes of an object for which implicit objectRegistry management is supported,...
void transfer(List< T > &list)
void resize(const label len)
Adjust allocated size of list.
MeshedSurface subsetMesh(const UList< bool > &include, labelList &pointMap, labelList &faceMap) const
Return a new surface subsetted on the selected faces.
const surfZoneList & surfZones() const
Const access to the surface zones.
void transfer(pointField &pointLst, List< Face > &faceLst)
Transfer the components.
label size() const
The surface size is the number of faces.
virtual void clear()
Clear all storage.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
const word & constant() const
Return constant name.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
void size(const label n)
Older name for setAddressableSize.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
void set(const bitSet &bitset)
Set specified bits from another bitset.
bool none() const
True if no bits in this bitset are set.
bool all() const
True if all bits in this bitset are set or if the set is empty.
const point & min() const
Minimum describing the bounding box.
const point & max() const
Maximum describing the bounding box.
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
bool intersect(const boundBox &bb)
Intersection (union) with the second box.
vector span() const
The bounding box span (from minimum to maximum)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Abstract base class for volume field interpolation.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
scalar print()
Print to screen.
Mesh consisting of general polyhedral cells.
virtual const faceList & faces() const
Return raw faces.
const boundBox & bounds() const
Return mesh bounding box.
const vectorField & faceCentres() const
const vectorField & cellCentres() const
A sampledSurface from a meshed surface. It samples on the points/faces of the meshed surface.
virtual bool expire()
Mark the surface as needing an update.
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool update()
Update the surface as required.
samplingSource
Types of sampling regions.
An abstract class for surfaces with sampling.
bool isPointData() const noexcept
Using interpolation to surface points.
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
bool interpolate() const noexcept
Same as isPointData()
A surface zone on a MeshedSurface.
labelRange range() const
The start/size range of this zone.
A class for managing temporary objects.
Standard boundBox with extra functionality for use in octree.
A List of wordRe with additional matching capabilities.
static wordRes uniq(const UList< wordRe > &input)
Return a wordRes with duplicate entries filtered out.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
OBJstream os(runTime.globalPath()/outputName)
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
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))
#define WarningInFunction
Report a warning using Foam::Warning.
labelList findMatching(const StringListType &input, const wordRes::filter &pred, AccessOp aop=identityOp())
Return ids for items with matching names.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
vector point
Point is a vector.
static IOobject selectReadIO(const word &name, const Time &runTime)
Ostream & endl(Ostream &os)
Add newline and flush stream.
MeshedSurface< face > meshedSurface
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
PointHit< point > pointHit
A PointIndexHit for 3D points.
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
UIndirectList< label > labelUIndList
UIndirectList of labels.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.
Extract name (as a word) from an object, typically using its name() method.