57template<
class ParcelType>
60template<
class ParcelType>
71template<
class ParcelType>
85 public ParcelType::constantProperties
128 inline scalar
T0()
const;
131 inline scalar
TMin()
const;
134 inline scalar
TMax()
const;
141 inline scalar
Cp0()
const;
149 inline scalar
f0()
const;
155 public ParcelType::trackingData
196 typedef typename ParcelType::trackingData::trackPart
trackPart;
201 template <
class TrackCloudType>
204 const TrackCloudType&
cloud,
205 trackPart part = ParcelType::trackingData::tpLinearTrack
234 inline scalar
Tc()
const;
240 inline scalar
Cpc()
const;
243 inline scalar&
Cpc();
263 template<
class TrackCloudType>
266 TrackCloudType&
cloud,
304 const label tetFacei,
323 const label tetFacei,
326 const scalar nParticle0,
328 const scalar dTarget0,
331 const vector& angularMomentum0,
342 bool newFormat =
true
366 const polyMesh& mesh_;
390 inline scalar
T()
const;
393 inline scalar
Cp()
const;
396 inline scalar
hs()
const;
411 template<
class TrackCloudType>
415 template<
class TrackCloudType>
418 TrackCloudType&
cloud,
424 template<
class TrackCloudType>
427 TrackCloudType&
cloud,
438 template<
class TrackCloudType>
441 TrackCloudType&
cloud,
450 template<
class CloudType>
454 template<
class CloudType>
463 const bool namesOnly =
false
467 template<
class CloudType>
471 template<
class CloudType>
477 friend Ostream& operator<< <ParcelType>
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
scalar Re(const trackingData &td) const
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Class to hold thermo particle constant properties.
constantProperties()
Null constructor.
void setTMax(const scalar TMax)
Set the maximum temperature [K].
scalar f0() const
Return const access to the particle scattering factor [].
scalar TMin() const
Return const access to minimum temperature [K].
scalar epsilon0() const
Return const access to the particle emissivity [].
scalar Cp0() const
Return const access to the particle specific heat capacity.
scalar TMax() const
Return const access to maximum temperature [K].
scalar T0() const
Return const access to the particle initial temperature [K].
iNew(const polyMesh &mesh)
autoPtr< ThermoParcel< ParcelType > > operator()(Istream &is) const
const interpolation< scalar > & GInterp() const
Return const access to the interpolator for continuous.
const volScalarField & kappa() const
Return access to the locally stored carrier kappa field.
const interpolation< scalar > & CpInterp() const
Return const access to the interpolator for continuous.
scalar Cpc() const
Return the continuous phase specific heat capacity.
const interpolation< scalar > & TInterp() const
Return const access to the interpolator for continuous.
const interpolation< scalar > & kappaInterp() const
Return const access to the interpolator for continuous.
ParcelType::trackingData::trackPart trackPart
const volScalarField & Cp() const
Return access to the locally stored carrier Cp field.
trackingData(const TrackCloudType &cloud, trackPart part=ParcelType::trackingData::tpLinearTrack)
Construct from components.
scalar Tc() const
Return the continuous phase temperature.
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
scalar Cp() const
Return const access to specific heat capacity.
scalar calcHeatTransfer(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar kappa, const scalar NCpW, const scalar Sh, scalar &dhsTrans, scalar &Sph)
Calculate new particle temperature.
scalar T() const
Return const access to temperature.
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Return a (basic particle) clone.
AddToPropertyList(ParcelType, " T"+" Cp")
String representation of properties.
scalar hs() const
Return the parcel sensible enthalpy.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
scalar & T()
Return access to temperature.
ThermoParcel(const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
Construct from Istream.
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label typeId, const scalar nParticle0, const scalar d0, const scalar dTarget0, const vector &U0, const vector &f0, const vector &angularMomentum0, const vector &torque0, const constantProperties &constProps)
Construct from components.
static const std::size_t sizeofFields
ThermoParcel(const ThermoParcel &p, const polyMesh &mesh)
Construct as a copy.
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
ThermoParcel(const ThermoParcel &p)
Construct as a copy.
ThermoParcel(const polyMesh &mesh, const vector &position, const label celli)
Construct from a position and a cell, searching for the rest of the.
static void writeFields(const CloudType &c)
Write.
void calcSurfaceValues(TrackCloudType &cloud, trackingData &td, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const
Calculate surface thermo properties.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
scalar & Cp()
Return access to specific heat capacity.
static void readFields(CloudType &c)
TypeName("ThermoParcel")
Runtime type information.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A cloud is a registry collection of lagrangian particles.
Class for demand-driven dictionary entries.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base class for volume field interpolation.
Registry of regIOobjects.
vector position() const
Return current particle position.
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.
Mesh consisting of general polyhedral cells.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
PtrList< coordinateSystem > coordinates(solidRegions.size())
OBJstream os(runTime.globalPath()/outputName)
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
DSMCCloud< dsmcParcel > CloudType
GeometricField< scalar, fvPatchField, volMesh > volScalarField
scalarField Re(const UList< complex > &cmplx)
Extract real component.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.