48 void Foam::sampledSets::lineFace::calcSamples
54 const bool storeFaces,
55 const bool storeCells,
85 for (
label bHitj = bHiti + 1; bHitj < bHits.
size(); ++ bHitj)
88 (bHits[bHiti].hitPoint() + bHits[bHitj].hitPoint())/2;
91 const scalar midT =
mag(midP - start)/
mag(end - start);
124 forAll(procCandidateCells, proci)
126 forAll(procCandidateCells[proci], candidatei)
128 const label celli = procCandidateCells[proci][candidatei];
130 if (celli == -1)
continue;
132 const scalar t = procCandidateTs[proci][candidatei];
133 const point p = (1 - t)*start + t*end;
143 i == 0 ? endStart : startEnd,
147 halfSegmentPositions,
148 halfSegmentDistances,
170 particles.
move(particles, tdBwd, rootGreat);
176 ?
mag(particles.
first()->position() - start)
178 : i == 0 ? vGreat : -vGreat,
179 [i](
const scalar a,
const scalar
b)
181 return (i == 0) == (a <
b) ? a :
b;
195 samplingPositions.
append(halfSegmentPositions);
197 samplingCells.
append(halfSegmentCells);
198 samplingFaces.
append(halfSegmentFaces);
200 halfSegmentPositions.
clear();
201 halfSegmentDistances.
clear();
202 halfSegmentCells.
clear();
203 halfSegmentFaces.
clear();
209 particle trackBwd(mesh, p, celli), trackFwd(trackBwd);
210 trackBwd.trackToFace(start - p, 0);
212 if (trackBwd.onFace() && trackFwd.
onFace())
216 (trackBwd.position() + trackFwd.
position())/2
218 samplingSegments.
append(segmenti);
219 samplingCells.
append(celli);
227 forAll(procCandidateCells, procj)
229 forAll(procCandidateCells[procj], candidatej)
231 const label cellj = procCandidateCells[procj][candidatej];
233 if (cellj == -1)
continue;
235 const scalar t = procCandidateTs[procj][candidatej];
239 t > segmentT.
first() - rootSmall
240 && t < segmentT.
second() + rootSmall
243 procCandidateCells[procj][candidatej] = -1;
244 procCandidateTs[procj][candidatej] = NaN;
255 samplingDistances.
resize(samplingPositions.
size());
256 forAll(samplingPositions, i)
258 samplingDistances[i] =
mag(samplingPositions[i] - start);
265 void Foam::sampledSets::lineFace::calcSamples
291 void Foam::sampledSets::lineFace::genSamples()
308 samplingPositions.
shrink();
309 samplingDistances.
shrink();
310 samplingSegments.
shrink();
336 start_(dict.
lookup(
"start")),
Template class for intrusive linked lists.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Particle for generating line-type sampled sets.
#define forAll(list, i)
Loop across all elements in list.
void resize(const label)
Alter the addressed list size.
A list of keyword definitions, which are a keyword followed by any number of values (e...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
void size(const label)
Override size to be inconsistent with allocated storage.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
virtual ~lineFace()
Destructor.
label size() const
Return the number of particles in the cloud.
Macros for easy insertion into run-time selection tables.
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Holds list of sampling points which is filled at construction time. Various implementations of this b...
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
void clear()
Clear the list, i.e. set size to zero.
A class for handling words, derived from string.
void append(const T &)
Append an element at the end of the list.
A Cloud of sampledSet particles.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
List< label > labelList
A List of labels.
addToRunTimeSelectionTable(sampledSet, arcUniform, word)
void reverse(UList< T > &, const label n)
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
bool onFace() const
Is the particle on a face?
const Type & second() const
Return second.
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
dimensioned< scalar > mag(const dimensioned< Type > &)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
T * first()
Return the first entry.
defineTypeNameAndDebug(arcUniform, 0)
void clear()
Clear the addressed list, i.e. set the size to zero.
lineFace(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
const Type & first() const
Return first.
vector position() const
Return current particle position.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.