31 template<
class ParcelType>
32 Foam::string Foam::ReactingParcel<ParcelType>::propertyList_ =
35 template<
class ParcelType>
44 template<
class ParcelType>
52 ParcelType(mesh, is, readFields),
61 if (is.
format() == IOstream::ASCII)
67 is.
read(reinterpret_cast<char*>(&mass0_), sizeofFields_);
77 "ReactingParcel<ParcelType>::ReactingParcel" 87 template<
class ParcelType>
88 template<
class CloudType>
95 template<
class ParcelType>
96 template<
class CloudType,
class CompositionType>
100 const CompositionType& compModel
103 bool valid = c.
size();
122 const wordList& phaseTypes = compModel.phaseTypes();
125 if (compModel.nPhase() == 1)
127 stateLabels = compModel.stateLabels()[0];
135 p.
Y_.setSize(nPhases, 0.0);
145 "Y" + phaseTypes[j] + stateLabels[j],
161 template<
class ParcelType>
162 template<
class CloudType>
165 ParcelType::writeFields(c);
169 template<
class ParcelType>
170 template<
class CloudType,
class CompositionType>
174 const CompositionType& compModel
177 ParcelType::writeFields(c);
193 const wordList& phaseTypes = compModel.phaseTypes();
195 if (compModel.nPhase() == 1)
197 stateLabels = compModel.stateLabels()[0];
206 "Y" + phaseTypes[j] + stateLabels[j],
231 template<
class ParcelType>
238 if (os.format() == IOstream::ASCII)
240 os << static_cast<const ParcelType&>(
p)
241 << token::SPACE <<
p.mass0()
242 << token::SPACE <<
p.Y();
246 os << static_cast<const ParcelType&>(
p);
249 reinterpret_cast<const char*>(&
p.mass0_),
258 "Ostream& operator<<(Ostream&, const ReactingParcel<ParcelType>&)" #define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
void size(const label)
Override size to be inconsistent with allocated storage.
scalarField Y_
Mass fractions of mixture [].
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
virtual Istream & read(token &)=0
Return next token from stream.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
streamFormat format() const
Return current stream format.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Base cloud calls templated on particle type.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
PtrList< volScalarField > & Y
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Mesh consisting of general polyhedral cells.
Reacting parcel class with one/two-way coupling with the continuous phase.
scalar mass0_
Initial mass [kg].
void transfer(List< T > &)
Transfer contents of the argument List into this.
A class for handling character strings derived from std::string.
Templated base class for dsmc cloud.
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
A primitive field of type <T> with automated input and output.
const scalarField & Y() const
Return const access to mass fractions of mixture [].