Go to the documentation of this file.
34 template<
class ParcelType>
39 template<
class ParcelType>
48 template<
class ParcelType>
74 template<
class ParcelType>
75 template<
class CloudType>
82 template<
class ParcelType>
83 template<
class CloudType,
class CompositionType>
87 const CompositionType& compModel
90 const bool valid =
c.size();
96 c.fieldIOobject(
"mass0", IOobject::MUST_READ),
99 c.checkFieldIOobject(
c, mass0);
104 p.mass0() = mass0[i];
108 const label idSolid = compModel.idSolid();
118 p.F_.setSize(nF,
Zero);
121 for (
label i = 0; i < nF; i++)
165 template<
class ParcelType>
166 template<
class CloudType>
176 template<
class ParcelType>
177 template<
class CloudType,
class CompositionType>
181 const CompositionType& compModel
188 const label np =
c.size();
189 const bool valid = np;
207 for (
label i = 0; i < nF; i++)
227 const label idSolid = compModel.idSolid();
254 template<
class ParcelType>
263 ParcelType::writeProperties(os, filters, delim, namesOnly);
266 #define writeProp(Name, Value) \
267 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
276 template<
class ParcelType>
277 template<
class CloudType>
284 ParcelType::readObjects(
c, obr);
288 template<
class ParcelType>
289 template<
class CloudType>
296 ParcelType::writeObjects(
c, obr);
300 template<
class ParcelType>
301 template<
class CloudType,
class CompositionType>
305 const CompositionType& compModel,
316 <<
"Reading of objects is still a work-in-progress" <<
nl;
321 template<
class ParcelType>
322 template<
class CloudType,
class CompositionType>
326 const CompositionType& compModel,
334 const label np =
c.size();
346 for (
label i = 0; i < nF; i++)
348 const word fieldName =
"F" +
name(i);
349 auto&
F = cloud::createIOField<scalar>(fieldName, np, obr);
359 const label idSolid = compModel.idSolid();
364 auto&
Y = cloud::createIOField<scalar>(fieldName, np, obr);
379 template<
class ParcelType>
387 if (os.format() == IOstream::ASCII)
389 os << static_cast<const ParcelType&>(
p)
390 << token::SPACE <<
F;
394 os << static_cast<const ParcelType&>(
p);
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
A class for handling words, derived from Foam::string.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
A primitive field of type <T> with automated input and output.
static constexpr const zero Zero
Global zero.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
volVectorField F(fluid.F())
A class for handling character strings derived from std::string.
Mesh consisting of general polyhedral cells.
#define forAll(list, i)
Loop across all elements in list.
Registry of regIOobjects.
Reacting heterogeneous Parcel.
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)....
word name(const complex &c)
Return string representation of complex.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields)
Templated base class for dsmc cloud.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
ReactingHeterogeneousParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
PtrList< volScalarField > & Y
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
A List of wordRe with additional matching capabilities.
const dimensionedScalar c
Speed of light in a vacuum.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
const volScalarField & p0
const wordList solidNames(rp["solid"])
#define WarningInFunction
Report a warning using Foam::Warning.
static const std::size_t sizeofFields
Size in bytes of the fields.
#define writeProp(Name, Value)