Go to the documentation of this file.
50 scalar fraction = testParticle.
trackToFace(position()+dt*
U, td);
61 const point& position,
72 sampledScalars_.setSize(td.
vsInterp_.size());
75 sampledScalars_[scalarI].append
86 sampledVectors_.setSize(td.
vvInterp_.size());
89 sampledVectors_[vectorI].append
136 is >> lifeTime_ >> sampledPositions_ >> sampledScalars
142 sampledScalars_[i].transfer(sampledScalars[i]);
144 sampledVectors_.setSize(sampledVectors.size());
147 sampledVectors_[i].transfer(sampledVectors[i]);
154 "streamLineParticle::streamLineParticle"
155 "(const Cloud<streamLineParticle>&, Istream&, bool)"
166 lifeTime_(
p.lifeTime_),
167 sampledPositions_(
p.sampledPositions_),
168 sampledScalars_(
p.sampledScalars_)
177 const scalar trackTime
185 scalar tEnd = (1.0 - stepFraction())*trackTime;
186 scalar maxDt = mesh_.bounds().mag();
206 sampledPositions_.append(position());
214 scalar magU =
mag(
U);
234 dt = calcSubCycleDeltaT(td, dt,
U);
243 scalar fraction = trackToFace(position() + dt*
U, td);
247 stepFraction() = 1.0 - tEnd/trackTime;
249 if (tEnd <= ROOTVSMALL)
275 Pout<<
"streamLineParticle : Removing stagnant particle:"
277 <<
" sampled positions:" << sampledPositions_.size()
285 sampledPositions_.append(position());
286 interpolateFields(td, position(),
cell(),
face());
290 Pout<<
"streamLineParticle : Removing particle:"
292 <<
" sampled positions:" << sampledPositions_.size()
303 forAll(sampledScalars_, i)
310 forAll(sampledVectors_, i)
328 const scalar trackFraction,
428 c.checkFieldIOobject(
c, lifeTime);
434 c.checkFieldIOobject(
c, sampledPositions);
445 iter().lifeTime_ = lifeTime[i];
446 iter().sampledPositions_.transfer(sampledPositions[i]);
478 lifeTime[i] = iter().lifeTime_;
479 sampledPositions[i] = iter().sampledPositions_;
485 sampledPositions.
write();
494 os << static_cast<const particle&>(
p)
501 os.
check(
"Ostream& operator<<(Ostream&, const streamLineParticle&)");
List< scalar > scalarList
A List of scalars.
Class used to pass tracking data to the trackToFace function.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
A primitive field of type <T> with automated input and output.
Wedge front and back plane patch.
#define forAll(list, i)
Loop across all elements in list.
Particle class that samples fields as it passes through. Used in streamline calculation.
scalar trackToFace(const vector &endPosition, TrackData &td)
Track particle to a given position and returns 1.0 if the.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
static void writeFields(const Cloud< streamLineParticle > &)
Write.
virtual bool write() const
Write using setting from DB.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void hitSymmetryPatch(const symmetryPolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a.
List< DynamicList< scalarList > > & allScalars_
dimensioned< scalar > mag(const dimensioned< Type > &)
This function object reads fields from the time directories and adds them to the mesh database for fu...
DynamicList< vectorList > & allPositions_
Mesh consisting of general polyhedral cells.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
const scalar trackLength_
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A patch is a list of labels that address the faces in the global face list.
Symmetry patch for non-planar or multi-plane patches.
void hitWallPatch(const wallPolyPatch &, trackingData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
Neighbour processor patch.
bool hitPatch(const polyPatch &, trackingData &td, const label patchI, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a patch.
List< DynamicList< vectorList > > & allVectors_
Ostream & operator<<(Ostream &, const edgeMesh &)
const PtrList< interpolation< scalar > > & vsInterp_
virtual bool check(const char *operation) const
Check IOstream status for given operation.
errorManip< error > abort(error &err)
void setSize(const label)
Reset size of List.
bool switchProcessor
Flag to switch processor.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
List< vector > vectorList
A List of vectors.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
scalar calcSubCycleDeltaT(trackingData &td, const scalar dt, const vector &U) const
Estimate dt to cross from current face to next one in nSubCycle.
streamLineParticle(const polyMesh &c, const vector &position, const label cellI, const label lifeTime)
Construct from components.
prefixOSstream Pout(cout, "Pout")
void hitProcessorPatch(const processorPolyPatch &, trackingData &td)
const PtrList< interpolation< vector > > & vvInterp_
const dimensionedScalar c
Speed of light in a vacuum.
static void readFields(CloudType &c)
Read the fields associated with the owner cloud.
A face is a list of labels corresponding to mesh vertices.
void hitWedgePatch(const wedgePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a wedge.
void size(const label)
Override size to be inconsistent with allocated storage.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static void writeFields(const CloudType &c)
Write the fields associated with the owner cloud.
vector interpolateFields(const trackingData &, const point &, const label cellI, const label faceI)
Interpolate all quantities; return interpolated velocity.
A cell is defined as a list of faces with extra functionality.
void hitCyclicPatch(const cyclicPolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a cyclic.
bool move(trackingData &, const scalar trackTime)
Track all particles to their end point.
void hitSymmetryPlanePatch(const symmetryPlanePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a.
static void readFields(Cloud< streamLineParticle > &)
Read.