58 for (
int i = 0; i < pTraits<Type>::nComponents; ++i)
80 if (!filters.empty() && !filters.match(
name))
88 writePropertyName<Type>(
os,
name, delim);
102 const Field<Type>& values,
105 const wordRes& filters
108 if (!filters.empty() && !filters.match(
name))
123 writePropertyName<Type>(
os, tag, delim);
134 os << delim << values;
139template<
class TrackCloudType>
142 const bool readOnProc = c.size();
146 const bool haveFile = procIO.typeHeaderOk<IOField<label>>(
true);
148 IOField<label> origProcId(procIO, readOnProc && haveFile);
149 c.checkFieldIOobject(c, origProcId);
151 IOField<label> origId
154 readOnProc && haveFile
156 c.checkFieldIOobject(c, origId);
159 for (particle&
p : c)
161 p.origProc_ = origProcId[i];
162 p.origId_ = origId[i];
169template<
class TrackCloudType>
172 const label np = c.size();
173 const bool writeOnProc = c.size();
175 if (writeLagrangianCoordinates)
177 IOPosition<TrackCloudType> ioP(c);
178 ioP.write(writeOnProc);
180 else if (!writeLagrangianPositions)
183 <<
"Must select coordinates and/or positions" <<
nl
188 if (writeLagrangianPositions)
195 ioP.write(writeOnProc);
212 origProc[i] =
p.origProc_;
213 origId[i] =
p.origId_;
219 origId.write(writeOnProc);
223template<
class CloudType>
230 const label np = c.size();
231 const label newNp = (positionPtr ? positionPtr->size() : 0);
234 for (label i = newNp; i < np; ++i)
236 parcelType*
p = c.last();
238 c.deleteParticle(*
p);
243 const auto& position = *positionPtr;
249 for (label i = np; i < newNp; ++i)
251 c.addParticle(
new parcelType(c.pMesh(), position[i], -1));
257 p.origProc_ = origProcId[i];
258 p.origId_ = origId[i];
263 p.relocate(position[i]);
272template<
class CloudType>
275 const label np =
c.size();
284 origProc[i] =
p.origProc_;
285 origId[i] =
p.origId_;
286 position[i] =
p.position();
293template<
class TrackCloudType>
296 const vector& displacement,
297 TrackCloudType& cloud,
301 typename TrackCloudType::particleType&
p =
302 static_cast<typename TrackCloudType::particleType&
>(*this);
303 typename TrackCloudType::particleType::trackingData& ttd =
304 static_cast<typename TrackCloudType::particleType::trackingData&
>(
td);
306 if (!
p.hitPatch(cloud, ttd))
308 const polyPatch& patch = mesh_.boundaryMesh()[
p.patch()];
312 p.hitWedgePatch(cloud, ttd);
316 p.hitSymmetryPlanePatch(
cloud, ttd);
320 p.hitSymmetryPatch(
cloud, ttd);
324 p.hitCyclicPatch(
cloud, ttd);
328 p.hitCyclicACMIPatch(
cloud, ttd, displacement);
332 p.hitCyclicAMIPatch(
cloud, ttd, displacement);
336 p.hitProcessorPatch(
cloud, ttd);
340 p.hitWallPatch(
cloud, ttd);
344 td.keepParticle =
false;
350template<
class TrackCloudType>
353 const vector& displacement,
354 TrackCloudType&
cloud,
362 else if (onInternalFace())
366 else if (onBoundaryFace())
370 changeToMasterPatch();
377template<
class TrackCloudType>
381 const scalar fraction,
382 TrackCloudType& cloud,
392template<
class TrackCloudType>
399template<
class TrackCloudType>
403 <<
"Hitting a wedge patch should not be possible."
410template<
class TrackCloudType>
413 TrackCloudType&
cloud,
421template<
class TrackCloudType>
430template<
class TrackCloudType>
436 const label receiveFacei = receiveCpp.
whichFace(facei_);
440 celli_ = mesh_.faceOwner()[facei_];
442 tetPti_ = mesh_.faces()[tetFacei_].size() - 1 - tetPti_;
454 : receiveCpp.
forwardT()[receiveFacei]
456 transformProperties(
T);
471template<
class TrackCloudType>
476 const vector& displacement
481 const cyclicAMIPolyPatch& cpp =
482 static_cast<const cyclicAMIPolyPatch&
>(mesh_.boundaryMesh()[patch()]);
483 const cyclicAMIPolyPatch& receiveCpp = cpp.
neighbPatch();
484 const label sendFacei = cpp.
whichFace(facei_);
485 const label receiveFacei = cpp.pointFace(sendFacei, displacement,
pos);
487 if (receiveFacei < 0)
491 td.keepParticle =
false;
493 <<
"Particle lost across " << cyclicAMIPolyPatch::typeName
494 <<
" patches " << cpp.name() <<
" and " << receiveCpp.name()
495 <<
" at position " <<
pos <<
endl;
499 facei_ = tetFacei_ = receiveFacei + receiveCpp.start();
502 vector displacementT = displacement;
503 cpp.reverseTransformDirection(displacementT, sendFacei);
509 const vector dispDir = cpp.fraction()*displacementT;
510 stepFraction_ += cpp.fraction();
515 mesh_.faceOwner()[facei_],
517 "Particle crossed between " + cyclicAMIPolyPatch::typeName +
518 " patches " + cpp.name() +
" and " + receiveCpp.name() +
519 " to a location outside of the mesh."
527 if (!receiveCpp.parallel())
531 receiveCpp.forwardT().size() == 1
532 ? receiveCpp.forwardT()[0]
533 : receiveCpp.forwardT()[receiveFacei]
535 transformProperties(
T);
537 else if (receiveCpp.separated())
541 (receiveCpp.separation().size() == 1)
542 ? receiveCpp.separation()[0]
543 : receiveCpp.separation()[receiveFacei]
545 transformProperties(-
s);
570template<
class TrackCloudType>
573 TrackCloudType& cloud,
575 const vector& displacement
578 const cyclicACMIPolyPatch& cpp =
579 static_cast<const cyclicACMIPolyPatch&
>(mesh_.boundaryMesh()[patch()]);
581 const label localFacei = cpp.whichFace(facei_);
585 const scalar mask = cpp.mask()[localFacei];
586 bool couple = mask >= 1 - cpp.tolerance();
587 bool nonOverlap = mask <= cpp.tolerance();
592 if (!couple && !nonOverlap)
595 couple = cpp.pointFace(localFacei, displacement,
pos) >= 0;
596 nonOverlap = !couple;
601 hitCyclicAMIPatch(
cloud,
td, displacement);
607 tetFacei_ = facei_ = cpp.nonOverlapPatch().start() + localFacei;
615template<
class TrackCloudType>
620template<
class TrackCloudType>
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
A primitive field of type <T> with automated input and output.
Helper IO class to read and write particle coordinates (positions).
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
bool empty() const noexcept
True if List is empty (ie, size() is zero).
A cloud is a registry collection of lagrangian particles.
static const IOField< Type > & lookupIOField(const word &fieldName, const objectRegistry &obr)
Lookup an IOField within object registry.
static const IOField< point > * findIOPosition(const objectRegistry &obr)
Locate the "position" IOField within object registry.
static IOField< Type > & createIOField(const word &fieldName, const label nParticle, objectRegistry &obr)
Helper to construct IOField on a supplied object registry.
virtual bool separated() const
Are the planes separated.
virtual bool parallel() const
Are the cyclic planes parallel.
virtual const vectorField & separation() const
If the planes are separated the separation vector.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
const polyPatch & nonOverlapPatch() const
Return a const reference to the non-overlapping patch.
const scalarField & mask() const
Mask field where 1 = overlap(coupled), 0 = no-overlap.
static scalar tolerance()
Overlap tolerance.
Cyclic patch for Arbitrary Mesh Interface (AMI).
label pointFace(const label facei, const vector &n, point &p) const
Return face index on neighbour patch which shares point p following trajectory vector n.
scalar fraction() const
Particle fraction increase between AMI pathces.
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to nbr side.
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
label transformGlobalFace(const label facei) const
const cyclicPolyPatch & neighbPatch() const
Registry of regIOobjects.
A traits class, which is primarily used for primitives and vector-space.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
static void writeProperty(Ostream &os, const word &name, const Type &value, const bool nameOnly, const word &delim, const wordRes &filters=wordRes::null())
Write a named particle property to stream, optionally filtered based on its name.
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but also stops on internal faces.
void hitCyclicACMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a.
vector position() const
Return current particle position.
void hitFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
void hitBoundaryFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Dispatch function for boundary face interaction. Calls one of.
bool onFace() const noexcept
Is the particle on a face?
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a cyclicPatch.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
bool onBoundaryFace() const noexcept
Is the particle on a boundary face?
label patch() const
Return the index of patch that the particle is on.
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
void hitCyclicAMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a cyclicAMIPatch.
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a symmetryPatch.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
static bool writeLagrangianPositions
Write particle positions file (v1706 format and earlier) Default is true (disable in etc/controlDict)...
bool onInternalFace() const noexcept
Is the particle on an internal face?
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
label origId() const noexcept
Return the particle ID on the originating processor.
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
vector normal() const
The (unit) normal of the tri on tetFacei_ for the current tet.
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
static void writePropertyName(Ostream &os, const word &name, const word &delim)
Write the name representation to stream.
static bool writeLagrangianCoordinates
Write particle coordinates file (v1712 and later) Default is true.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
label origProc() const noexcept
Return the originating processor ID.
const word & name() const noexcept
The patch name.
A patch is a list of labels that address the faces in the global face list.
label whichFace(const label facei) const noexcept
Return label of face in patch from global face label.
label start() const noexcept
Return start label of this patch in the polyMesh face list.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A List of wordRe with additional matching capabilities.
static bool match(const UList< wordRe > &selectors, const std::string &text, bool literal=false)
Test for a match of any selectors against the text.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
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))
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
#define WarningInFunction
Report a warning using Foam::Warning.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar pos(const dimensionedScalar &ds)
DSMCCloud< dsmcParcel > CloudType
static const Identity< scalar > I
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
errorManip< error > abort(error &err)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.