37#ifndef Foam_DTRMParticle_H
38#define Foam_DTRMParticle_H
86 label transmissiveId_;
144 inline scalar&
Q(label celli);
175 const vector& targetPosition,
179 const label transmissiveId
188 const label tetFacei,
191 const vector& targetPosition,
194 const label transmissiveId
203 bool newFormat =
true
307 const bool namesOnly =
false
Base cloud calls templated on particle type.
iNew(const polyMesh &mesh)
autoPtr< DTRMParticle > operator()(Istream &is) const
Class used to pass tracking data to the trackToFace function.
const UPtrList< reflectionModel > & reflection() const
trackingData(Cloud< DTRMParticle > &spc, const interpolationCell< scalar > &aInterp, const interpolationCell< scalar > &eInterp, const interpolationCell< scalar > &EInterp, const interpolationCell< scalar > &TInterp, const interpolationCellPoint< vector > &nHatInterp, const labelField &, const UPtrList< reflectionModel > &, volScalarField &Q)
const labelField & relfectedCells() const
const interpolationCell< scalar > & aInterp() const
const interpolationCellPoint< vector > & nHatInterp() const
const interpolationCell< scalar > & EInterp() const
const interpolationCell< scalar > & TInterp() const
const interpolationCell< scalar > & eInterp() const
Discrete Transfer Radiation Model (DTRM) particle.
DTRMParticle(const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
Construct from Istream.
scalar dA() const noexcept
Return const access dA.
point & p1() noexcept
Return access to the target position.
virtual autoPtr< particle > clone() const
Return a clone.
AddToPropertyList(particle, " p0"+" p1"+" I0"+" I"+" dA"+" transmissiveId";)
String representation of properties.
scalar I() const noexcept
Return const access to the current intensity.
DTRMParticle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const vector &position, const vector &targetPosition, const scalar I, const scalar dA, const label transmissiveId)
Construct from components.
DTRMParticle(const polyMesh &mesh, const vector &position, const vector &targetPosition, const scalar I, const label cellI, const scalar dA, const label transmissiveId)
Construct from components, with searching for tetFace and.
const point & p0() const noexcept
Return const access to the initial position.
friend Ostream & operator<<(Ostream &os, const DTRMParticle &p)
scalar I0() const noexcept
Return const access to the initial intensity.
DTRMParticle(const DTRMParticle &p)
Construct as copy.
const point & p1() const noexcept
Return const access to the target position.
scalar & I0() noexcept
Return access to the initial intensity.
void hitProcessorPatch(Cloud< DTRMParticle > &, trackingData &td)
Overridable function to handle the particle hitting a processorPatch.
static const std::size_t sizeofFields_
Size in bytes of the fields.
bool hitPatch(Cloud< DTRMParticle > &, trackingData &td)
bool move(Cloud< DTRMParticle > &, trackingData &, const scalar)
Move.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
void hitWallPatch(Cloud< DTRMParticle > &, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
scalar & I() noexcept
Return access to the current intensity.
scalar & dA() noexcept
Return access to dA.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Uses the cell value for any location within the cell.
vector position() const
Return current particle position.
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
const polyMesh & mesh() const noexcept
Return the mesh database.
const barycentric & coordinates() const noexcept
Return current particle coordinates.
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Mesh consisting of general polyhedral cells.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for radiation modelling.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Field< label > labelField
Specialisation of Field<T> for label.
vector point
Point is a vector.
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
Forwards and collection of common volume field types.