34template<
class ParcelType>
39template<
class ParcelType>
48template<
class ParcelType>
74template<
class ParcelType>
75template<
class CloudType>
78 ParcelType::readFields(c);
82template<
class ParcelType>
83template<
class CloudType,
class CompositionType>
87 const CompositionType& compModel
90 const bool readOnProc = c.size();
92 ParcelType::readFields(c);
99 c.checkFieldIOobject(c, mass0);
104 p.mass0() = mass0[i];
108 const label idSolid = compModel.idSolid();
115 for (ReactingHeterogeneousParcel<ParcelType>&
p : c)
118 p.F_.setSize(nF,
Zero);
121 for (label i = 0; i < nF; i++)
165template<
class ParcelType>
166template<
class CloudType>
172 ParcelType::writeFields(c);
176template<
class ParcelType>
177template<
class CloudType,
class CompositionType>
181 const CompositionType& compModel
188 const label np = c.size();
189 const bool writeOnProc = c.size();
205 mass0.write(writeOnProc);
207 for (label i = 0; i < nF; i++)
227 const label idSolid = compModel.idSolid();
249 Y.write(writeOnProc);
254template<
class ParcelType>
258 const wordRes& filters,
263 ParcelType::writeProperties(
os, filters, delim, namesOnly);
266 #define writeProp(Name, Value) \
267 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
276template<
class ParcelType>
277template<
class CloudType>
284 ParcelType::readObjects(c, obr);
288template<
class ParcelType>
289template<
class CloudType>
296 ParcelType::writeObjects(c, obr);
300template<
class ParcelType>
301template<
class CloudType,
class CompositionType>
305 const CompositionType& compModel,
316 <<
"Reading of objects is still a work-in-progress" <<
nl;
321template<
class ParcelType>
322template<
class CloudType,
class CompositionType>
326 const CompositionType& compModel,
334 const label np = c.size();
335 const bool writeOnProc = c.size();
347 for (label i = 0; i < nF; i++)
349 const word fieldName =
"F" +
name(i);
353 for (
const ReactingHeterogeneousParcel<ParcelType>&
p0 : c)
360 const label idSolid = compModel.idSolid();
380template<
class ParcelType>
381Foam::Ostream& Foam::operator<<
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
A primitive field of type <T> with automated input and output.
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
streamFormat format() const noexcept
Get the current stream format.
@ ASCII
"ascii" (normal default)
virtual bool check(const char *operation) const
Check IOstream status for given operation.
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,...
Reacting heterogeneous Parcel.
label canCombust_
Flag to identify if the particle can devolatilise and combust.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
ReactingHeterogeneousParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
static const std::size_t sizeofFields
Size in bytes of the fields.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
scalarField F_
Progress variables for reactions.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
const scalarField & F() const
Return const access to F.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
static IOField< Type > & createIOField(const word &fieldName, const label nParticle, objectRegistry &obr)
Helper to construct IOField on a supplied object registry.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A class for handling character strings derived from std::string.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
const volScalarField & p0
const wordList solidNames(rp["solid"])
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define WarningInFunction
Report a warning using Foam::Warning.
volVectorField F(fluid.F())
const dimensionedScalar c
Speed of light in a vacuum.
List< word > wordList
List of word.
DSMCCloud< dsmcParcel > CloudType
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
static constexpr const zero Zero
Global zero (0).
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.