ReactingParcelIO.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "ReactingParcel.H"
30#include "IOstreams.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34template<class ParcelType>
37
38
39template<class ParcelType>
41(
42 sizeof(scalar)
43);
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
48template<class ParcelType>
50(
51 const polyMesh& mesh,
52 Istream& is,
53 bool readFields,
54 bool newFormat
55)
56:
57 ParcelType(mesh, is, readFields, newFormat),
58 mass0_(0.0),
59 Y_(0)
60{
61 if (readFields)
62 {
64
65 if (is.format() == IOstream::ASCII)
66 {
67 is >> mass0_;
68 }
69 else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
70 {
71 // Non-native label or scalar size
72
73 is.beginRawRead();
74
75 readRawScalar(is, &mass0_);
76
77 is.endRawRead();
78 }
79 else
80 {
81 is.read(reinterpret_cast<char*>(&mass0_), sizeofFields);
82 }
83
84 is >> Ymix;
85
86 Y_.transfer(Ymix);
87 }
88
90}
91
92
93template<class ParcelType>
94template<class CloudType>
96{
98}
99
100
101template<class ParcelType>
102template<class CloudType, class CompositionType>
104(
105 CloudType& c,
106 const CompositionType& compModel
107)
108{
109 bool valid = c.size();
110
112
113 IOField<scalar> mass0
114 (
115 c.fieldIOobject("mass0", IOobject::MUST_READ),
116 valid
117 );
118 c.checkFieldIOobject(c, mass0);
119
120 label i = 0;
122 {
123 p.mass0_ = mass0[i];
124
125 ++i;
126 }
127
128 // Get names and sizes for each Y...
129 const wordList& phaseTypes = compModel.phaseTypes();
130 const label nPhases = phaseTypes.size();
131 wordList stateLabels(nPhases, "");
132 if (compModel.nPhase() == 1)
133 {
134 stateLabels = compModel.stateLabels()[0];
135 }
136
137
138 // Set storage for each Y... for each parcel
140 {
141 p.Y_.setSize(nPhases, 0.0);
142 }
143
144 // Populate Y for each parcel
145 forAll(phaseTypes, j)
146 {
148 (
149 c.fieldIOobject
150 (
151 "Y" + phaseTypes[j] + stateLabels[j],
153 ),
154 valid
155 );
156
157 label i = 0;
159 {
160 p.Y_[j] = Y[i];
161
162 ++i;
163 }
164 }
165}
166
167
168template<class ParcelType>
169template<class CloudType>
171{
173}
174
175
176template<class ParcelType>
177template<class CloudType, class CompositionType>
179(
180 const CloudType& c,
181 const CompositionType& compModel
182)
183{
185
186 const label np = c.size();
187
188 {
189 IOField<scalar> mass0(c.fieldIOobject("mass0", IOobject::NO_READ), np);
190
191 label i = 0;
192 for (const ReactingParcel<ParcelType>& p : c)
193 {
194 mass0[i] = p.mass0_;
195
196 ++i;
197 }
198 mass0.write(np > 0);
199
200 // Write the composition fractions
201 const wordList& phaseTypes = compModel.phaseTypes();
202 wordList stateLabels(phaseTypes.size(), "");
203 if (compModel.nPhase() == 1)
204 {
205 stateLabels = compModel.stateLabels()[0];
206 }
207
208 forAll(phaseTypes, j)
209 {
211 (
212 c.fieldIOobject
213 (
214 "Y" + phaseTypes[j] + stateLabels[j],
216 ),
217 np
218 );
219
220 label i = 0;
221 for (const ReactingParcel<ParcelType>& p : c)
222 {
223 Y[i] = p.Y()[j];
224
225 ++i;
226 }
227
228 Y.write(np > 0);
229 }
230 }
231}
232
233
234template<class ParcelType>
236(
237 Ostream& os,
238 const wordRes& filters,
239 const word& delim,
240 const bool namesOnly
241) const
242{
243 ParcelType::writeProperties(os, filters, delim, namesOnly);
244
245 #undef writeProp
246 #define writeProp(Name, Value) \
247 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
248
249 writeProp("mass0", mass0_);
250 writeProp("Y", Y_);
251
252 #undef writeProp
253}
254
255
256template<class ParcelType>
257template<class CloudType>
259(
260 CloudType& c,
261 const objectRegistry& obr
262)
263{
265}
266
267
268template<class ParcelType>
269template<class CloudType>
271(
272 const CloudType& c,
273 objectRegistry& obr
274)
275{
277}
278
279
280template<class ParcelType>
281template<class CloudType, class CompositionType>
283(
284 CloudType& c,
285 const CompositionType& compModel,
286 const objectRegistry& obr
287)
288{
290
291 if (!c.size()) return;
292
293
294 auto& mass0 = cloud::lookupIOField<scalar>("mass0", obr);
295
296 label i = 0;
298 {
299 p.mass0_ = mass0[i];
300
301 ++i;
302 }
303
304 // The composition fractions
305 const wordList& phaseTypes = compModel.phaseTypes();
306 wordList stateLabels(phaseTypes.size(), "");
307 if (compModel.nPhase() == 1)
308 {
309 stateLabels = compModel.stateLabels()[0];
310 }
311
312 forAll(phaseTypes, j)
313 {
314 const word fieldName = "Y" + phaseTypes[j] + stateLabels[j];
315 auto& Y = cloud::lookupIOField<scalar>(fieldName, obr);
316
317 label i = 0;
319 {
320 p.Y()[j] = Y[i];
321
322 ++i;
323 }
324 }
325}
326
327
328template<class ParcelType>
329template<class CloudType, class CompositionType>
331(
332 const CloudType& c,
333 const CompositionType& compModel,
334 objectRegistry& obr
335)
336{
338
339 const label np = c.size();
340
341 if (np > 0)
342 {
343 auto& mass0 = cloud::createIOField<scalar>("mass0", np, obr);
344
345 label i = 0;
346 for (const ReactingParcel<ParcelType>& p : c)
347 {
348 mass0[i] = p.mass0_;
349
350 ++i;
351 }
352
353 // Write the composition fractions
354 const wordList& phaseTypes = compModel.phaseTypes();
355 wordList stateLabels(phaseTypes.size(), "");
356 if (compModel.nPhase() == 1)
357 {
358 stateLabels = compModel.stateLabels()[0];
359 }
360
361 forAll(phaseTypes, j)
362 {
363 const word fieldName = "Y" + phaseTypes[j] + stateLabels[j];
364 auto& Y = cloud::createIOField<scalar>(fieldName, np, obr);
365
366 label i = 0;
367 for (const ReactingParcel<ParcelType>& p : c)
368 {
369 Y[i] = p.Y()[j];
370
371 ++i;
372 }
373 }
374 }
375}
376
377
378// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
379
380template<class ParcelType>
381Foam::Ostream& Foam::operator<<
382(
383 Ostream& os,
385)
386{
387 if (os.format() == IOstream::ASCII)
388 {
389 os << static_cast<const ParcelType&>(p)
390 << token::SPACE << p.mass0()
391 << token::SPACE << p.Y();
392 }
393 else
394 {
395 os << static_cast<const ParcelType&>(p);
396 os.write
397 (
398 reinterpret_cast<const char*>(&p.mass0_),
400 );
401 os << p.Y();
402 }
403
405 return os;
406}
407
408
409// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
streamFormat format() const noexcept
Get the current stream format.
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const noexcept
Definition: IOstream.H:300
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const noexcept
Definition: IOstream.H:309
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
virtual bool endRawRead()=0
End of low-level raw binary read.
virtual bool beginRawRead()=0
Start of low-level raw binary read.
virtual Istream & read(token &)=0
Return next token from stream.
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Reacting parcel class with one/two-way coupling with the continuous phase.
static const std::size_t sizeofFields
Size in bytes of the fields.
scalarField Y_
Mass fractions of mixture [].
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
scalar mass0_
Initial mass [kg].
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
Allows specification of different writing frequency of objects registered to the database.
Definition: writeObjects.H:142
static void readObjects(Cloud< injectedParticle > &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
Registry of regIOobjects.
static string propertyList()
Definition: particle.H:372
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for handling character strings derived from std::string.
Definition: string.H:79
@ SPACE
Space [isspace].
Definition: token.H:125
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define FUNCTION_NAME
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333