Go to the documentation of this file.
34 const point& position,
57 !sampledPositions_.size()
62 sampledPositions_.append(position);
65 sampledScalars_.setSize(td.
vsInterp_.size());
68 sampledScalars_[scalarI].append
80 sampledVectors_.setSize(td.
vvInterp_.size());
90 positionU = td.
vvInterp_[vectorI].interpolate
97 sampledVectors_[vectorI].append(positionU);
110 vector U = interpolateFields(td, position(),
cell(), tetFace());
117 scalar magU =
mag(
U);
138 const label tetFaceI,
140 const label meshEdgeStart,
141 const label diagEdge,
174 >> sampledPositions_ >> sampledScalars >> sampledVectors;
179 sampledScalars_[i].transfer(sampledScalars[i]);
181 sampledVectors_.setSize(sampledVectors.size());
184 sampledVectors_[i].transfer(sampledVectors[i]);
191 "wallBoundedStreamLineParticle::wallBoundedStreamLineParticle"
192 "(const Cloud<wallBoundedStreamLineParticle>&, Istream&, bool)"
203 lifeTime_(
p.lifeTime_),
204 sampledPositions_(
p.sampledPositions_),
205 sampledScalars_(
p.sampledScalars_),
206 sampledVectors_(
p.sampledVectors_)
215 const scalar trackTime
230 scalar tEnd = (1.0 - stepFraction())*trackTime;
231 scalar maxDt = mesh_.bounds().mag();
263 scalar fraction = trackToEdge(td, position() + dt*
U);
267 stepFraction() = 1.0 - tEnd/trackTime;
270 if (tEnd <= ROOTVSMALL)
294 Pout<<
"wallBoundedStreamLineParticle :"
295 <<
" Removing stagnant particle:"
297 <<
" sampled positions:" << sampledPositions_.size()
309 Pout<<
"wallBoundedStreamLineParticle : Removing particle:"
311 <<
" sampled positions:" << sampledPositions_.size()
324 forAll(sampledScalars_, i)
331 forAll(sampledVectors_, i)
358 c.fieldIOobject(
"lifeTime", IOobject::MUST_READ)
360 c.checkFieldIOobject(
c, lifeTime);
364 c.fieldIOobject(
"sampledPositions", IOobject::MUST_READ)
366 c.checkFieldIOobject(
c, sampledPositions);
371 iter().lifeTime_ = lifeTime[i];
372 iter().sampledPositions_.transfer(sampledPositions[i]);
383 wallBoundedParticle::writeFields(
c);
389 c.fieldIOobject(
"lifeTime", IOobject::NO_READ),
394 c.fieldIOobject(
"sampledPositions", IOobject::NO_READ),
401 lifeTime[i] = iter().lifeTime_;
402 sampledPositions[i] = iter().sampledPositions_;
407 sampledPositions.
write();
419 os << static_cast<const wallBoundedParticle&>(
p)
420 << token::SPACE <<
p.lifeTime_
421 << token::SPACE <<
p.sampledPositions_
422 << token::SPACE <<
p.sampledScalars_
423 << token::SPACE <<
p.sampledVectors_;
428 "Ostream& operator<<(Ostream&, const wallBoundedStreamLineParticle&)"
List< scalar > scalarList
A List of scalars.
bool move(trackingData &, const scalar trackTime)
Track all particles to their end point.
#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.
#define forAll(list, i)
Loop across all elements in list.
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
List< DynamicList< scalarList > > & allScalars_
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
virtual bool write() const
Write using setting from DB.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
This function object reads fields from the time directories and adds them to the mesh database for fu...
Mesh consisting of general polyhedral cells.
vector sample(trackingData &td)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
List< DynamicList< vectorList > > & allVectors_
DynamicList< vectorList > & allPositions_
Class used to pass tracking data to the trackToEdge function.
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)....
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
Particle class that samples fields as it passes through. Used in streamline calculation.
const PtrList< interpolation< vector > > & vvInterp_
static void writeFields(const Cloud< wallBoundedStreamLineParticle > &)
Write.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
vector interpolateFields(const trackingData &td, const point &position, const label cellI, const label faceI)
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.
static void readFields(Cloud< wallBoundedStreamLineParticle > &)
Read.
List< vector > vectorList
A List of vectors.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
prefixOSstream Pout(cout, "Pout")
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
const scalar trackLength_
wallBoundedStreamLineParticle(const polyMesh &c, const vector &position, const label cellI, const label tetFaceI, const label tetPtI, const label meshEdgeStart, const label diagEdge, const label lifeTime)
Construct from components.
const PtrList< interpolation< scalar > > & vsInterp_
const dimensionedScalar c
Speed of light in a vacuum.
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,...
A cell is defined as a list of faces with extra functionality.
dimensioned< scalar > magSqr(const dimensioned< Type > &)