Go to the documentation of this file.
48 void Foam::triSurfaceMeshPointSet::calcSamples
50 DynamicList<point>& samplingPts,
51 DynamicList<label>& samplingCells,
52 DynamicList<label>& samplingFaces,
53 DynamicList<label>& samplingSegments,
54 DynamicList<scalar>& samplingCurveDist
57 forAll(sampleCoords_, sampleI)
63 samplingPts.append(sampleCoords_[sampleI]);
64 samplingCells.append(celli);
65 samplingFaces.append(-1);
66 samplingSegments.append(0);
67 samplingCurveDist.append(1.0 * sampleI);
73 void Foam::triSurfaceMeshPointSet::genSamples()
76 DynamicList<point> samplingPts;
77 DynamicList<label> samplingCells;
78 DynamicList<label> samplingFaces;
79 DynamicList<label> samplingSegments;
80 DynamicList<scalar> samplingCurveDist;
92 samplingCells.shrink();
93 samplingFaces.shrink();
94 samplingSegments.shrink();
95 samplingCurveDist.shrink();
100 std::move(samplingPts),
101 std::move(samplingCells),
102 std::move(samplingFaces),
103 std::move(samplingSegments),
104 std::move(samplingCurveDist)
129 const auto* surfPtr =
136 sampleCoords_ = surfPtr->
points();
int debug
Static debugging option.
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A class for handling words, derived from Foam::string.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
static constexpr const zero Zero
Global zero (0)
const meshSearch & searchEngine() const
triSurfaceMeshPointSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
IOoject and searching on triSurface.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
virtual point getRefPoint(const List< point > &pts) const
Get reference point.
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Mesh consisting of general polyhedral cells.
#define forAll(list, i)
Loop across all elements in list.
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Triangulated surface description with patch information.
messageStream Info
Information stream (uses stdout - output is on the master only)
word name(const complex &c)
Return string representation of complex.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const Field< point_type > & points() const
Return reference to global points.
Macros for easy insertion into run-time selection tables.
virtual tmp< pointField > points() const
Get the points that define the surface.
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
const Time & time() const
Return the top-level database.
const word & constant() const
Return constant name.
defineTypeNameAndDebug(combustionModel, 0)