Loading...
Searching...
No Matches
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-2022 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() == IOstreamOption::ASCII)
66 {
67 is >> mass0_;
68 }
69 else
70 {
71 // No non-native streaming
73
74 is.read(reinterpret_cast<char*>(&mass0_), sizeofFields);
75 }
76
77 is >> Ymix;
78
79 Y_.transfer(Ymix);
80 }
83}
84
85
86template<class ParcelType>
87template<class CloudType>
90 ParcelType::readFields(c);
91}
92
93
94template<class ParcelType>
95template<class CloudType, class CompositionType>
97(
98 CloudType& c,
99 const CompositionType& compModel
100)
101{
102 const bool readOnProc = c.size();
103
104 ParcelType::readFields(c);
105
106 IOField<scalar> mass0
107 (
108 c.newIOobject("mass0", IOobject::MUST_READ),
109 readOnProc
110 );
111 c.checkFieldIOobject(c, mass0);
112
113 label i = 0;
115 {
116 p.mass0_ = mass0[i];
117
118 ++i;
119 }
120
121 // Get names and sizes for each Y...
122 const wordList& phaseTypes = compModel.phaseTypes();
123 const label nPhases = phaseTypes.size();
124 wordList stateLabels(nPhases, "");
125 if (compModel.nPhase() == 1)
126 {
127 stateLabels = compModel.stateLabels()[0];
128 }
129
130
131 // Set storage for each Y... for each parcel
133 {
134 p.Y_.setSize(nPhases, 0.0);
135 }
136
137 // Populate Y for each parcel
138 forAll(phaseTypes, j)
139 {
141 (
142 c.newIOobject
143 (
144 "Y" + phaseTypes[j] + stateLabels[j],
146 ),
147 readOnProc
148 );
149
150 label i = 0;
152 {
153 p.Y_[j] = Y[i];
154
155 ++i;
157 }
158}
159
160
161template<class ParcelType>
162template<class CloudType>
165 ParcelType::writeFields(c);
166}
167
168
169template<class ParcelType>
170template<class CloudType, class CompositionType>
172(
173 const CloudType& c,
174 const CompositionType& compModel
175)
176{
177 ParcelType::writeFields(c);
178
179 const label np = c.size();
180 const bool writeOnProc = c.size();
181
182 {
183 IOField<scalar> mass0(c.newIOobject("mass0", IOobject::NO_READ), np);
184
185 label i = 0;
186 for (const ReactingParcel<ParcelType>& p : c)
187 {
188 mass0[i] = p.mass0_;
189
190 ++i;
191 }
192 mass0.write(writeOnProc);
193
194 // Write the composition fractions
195 const wordList& phaseTypes = compModel.phaseTypes();
196 wordList stateLabels(phaseTypes.size(), "");
197 if (compModel.nPhase() == 1)
198 {
199 stateLabels = compModel.stateLabels()[0];
200 }
201
202 forAll(phaseTypes, j)
203 {
205 (
206 c.newIOobject
207 (
208 "Y" + phaseTypes[j] + stateLabels[j],
210 ),
211 np
212 );
213
214 label i = 0;
215 for (const ReactingParcel<ParcelType>& p : c)
216 {
217 Y[i] = p.Y()[j];
218
219 ++i;
220 }
221
222 Y.write(writeOnProc);
223 }
224 }
225}
226
227
228template<class ParcelType>
230(
231 Ostream& os,
232 const wordRes& filters,
233 const word& delim,
234 const bool namesOnly
235) const
236{
237 ParcelType::writeProperties(os, filters, delim, namesOnly);
238
239 #undef writeProp
240 #define writeProp(Name, Value) \
241 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
242
243 writeProp("mass0", mass0_);
244 writeProp("Y", Y_);
246 #undef writeProp
247}
248
249
250template<class ParcelType>
251template<class CloudType>
253(
254 CloudType& c,
255 const objectRegistry& obr
256)
258 ParcelType::readObjects(c, obr);
259}
260
261
262template<class ParcelType>
263template<class CloudType>
265(
266 const CloudType& c,
267 objectRegistry& obr
268)
270 ParcelType::writeObjects(c, obr);
271}
272
273
274template<class ParcelType>
275template<class CloudType, class CompositionType>
277(
278 CloudType& c,
279 const CompositionType& compModel,
280 const objectRegistry& obr
281)
282{
283 ParcelType::readObjects(c, obr);
284
285 if (!c.size()) return;
286
287
288 auto& mass0 = cloud::lookupIOField<scalar>("mass0", obr);
289
290 label i = 0;
292 {
293 p.mass0_ = mass0[i];
294
295 ++i;
296 }
297
298 // The composition fractions
299 const wordList& phaseTypes = compModel.phaseTypes();
300 wordList stateLabels(phaseTypes.size(), "");
301 if (compModel.nPhase() == 1)
302 {
303 stateLabels = compModel.stateLabels()[0];
304 }
305
306 forAll(phaseTypes, j)
307 {
308 const word fieldName = "Y" + phaseTypes[j] + stateLabels[j];
309 auto& Y = cloud::lookupIOField<scalar>(fieldName, obr);
310
311 label i = 0;
313 {
314 p.Y()[j] = Y[i];
315
316 ++i;
318 }
319}
320
321
322template<class ParcelType>
323template<class CloudType, class CompositionType>
325(
326 const CloudType& c,
327 const CompositionType& compModel,
328 objectRegistry& obr
329)
330{
331 ParcelType::writeObjects(c, obr);
332
333 const label np = c.size();
334 const bool writeOnProc = c.size();
335
336 if (writeOnProc)
337 {
338 auto& mass0 = cloud::createIOField<scalar>("mass0", np, obr);
339
340 label i = 0;
341 for (const ReactingParcel<ParcelType>& p : c)
342 {
343 mass0[i] = p.mass0_;
344
345 ++i;
346 }
347
348 // Write the composition fractions
349 const wordList& phaseTypes = compModel.phaseTypes();
350 wordList stateLabels(phaseTypes.size(), "");
351 if (compModel.nPhase() == 1)
352 {
353 stateLabels = compModel.stateLabels()[0];
354 }
355
356 forAll(phaseTypes, j)
357 {
358 const word fieldName = "Y" + phaseTypes[j] + stateLabels[j];
359 auto& Y = cloud::createIOField<scalar>(fieldName, np, obr);
360
361 label i = 0;
362 for (const ReactingParcel<ParcelType>& p : c)
363 {
364 Y[i] = p.Y()[j];
365
366 ++i;
367 }
368 }
369 }
370}
371
372
373// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
374
375template<class ParcelType>
376Foam::Ostream& Foam::operator<<
377(
378 Ostream& os,
380)
381{
383 {
385 << token::SPACE << p.mass0()
386 << token::SPACE << p.Y();
387 }
388 else
389 {
391 os.write
392 (
393 reinterpret_cast<const char*>(&p.mass0_),
395 );
396 os << p.Y();
397 }
398
400 return os;
401}
402
403
404// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
void setSize(const label n)
Alias for resize().
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ 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.
Definition IOstream.C:45
bool fatalCheckNativeSizes(const char *operation) const
Assert that the label/scalar byte-size associated with the stream are the native label/scalar sizes.
Definition IOstream.C:67
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
virtual Istream & read(token &)=0
Return next token from stream.
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
virtual Ostream & write(const char c) override
Write character.
Definition OBJstream.C:69
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Reacting parcel class with one/two-way coupling with the continuous phase.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
static const std::size_t sizeofFields
Size in bytes of the fields.
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 - composition supplied.
const scalarField & Y() const
Return const access to mass fractions of mixture [].
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() const
Return const access to initial mass [kg].
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:118
static const IOField< Type > & lookupIOField(const word &fieldName, const objectRegistry &obr)
Lookup an IOField within object registry.
Definition cloud.H:202
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.
Definition polyMesh.H:79
A class for handling character strings derived from std::string.
Definition string.H:76
@ SPACE
Space [isspace].
Definition token.H:144
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
PtrList< volScalarField > & Y
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define FUNCTION_NAME
const dimensionedScalar c
Speed of light in a vacuum.
List< word > wordList
List of word.
Definition fileName.H:60
DSMCCloud< dsmcParcel > CloudType
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299