ReactingMultiphaseParcelIO.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-2020 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
30#include "IOstreams.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34template<class ParcelType>
37
38
39template<class ParcelType>
41(
42 0
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 YGas_(0),
59 YLiquid_(0),
60 YSolid_(0),
61 canCombust_(0)
62{
63 if (readFields)
64 {
68
69 is >> Yg >> Yl >> Ys;
70
71 YGas_.transfer(Yg);
73 YSolid_.transfer(Ys);
74 }
75
77}
78
79
80template<class ParcelType>
81template<class CloudType>
83{
85}
86
87
88template<class ParcelType>
89template<class CloudType, class CompositionType>
91(
92 CloudType& c,
93 const CompositionType& compModel
94)
95{
96 bool valid = c.size();
97
98 ParcelType::readFields(c, compModel);
99
100 // Get names and sizes for each Y...
101 const label idGas = compModel.idGas();
102 const wordList& gasNames = compModel.componentNames(idGas);
103 const label idLiquid = compModel.idLiquid();
104 const wordList& liquidNames = compModel.componentNames(idLiquid);
105 const label idSolid = compModel.idSolid();
106 const wordList& solidNames = compModel.componentNames(idSolid);
107 const wordList& stateLabels = compModel.stateLabels();
108
109 // Set storage for each Y... for each parcel
111 {
112 p.YGas_.setSize(gasNames.size(), 0.0);
113 p.YLiquid_.setSize(liquidNames.size(), 0.0);
114 p.YSolid_.setSize(solidNames.size(), 0.0);
115 }
116
117 // Populate YGas for each parcel
118 forAll(gasNames, j)
119 {
120 IOField<scalar> YGas
121 (
122 c.fieldIOobject
123 (
124 "Y" + gasNames[j] + stateLabels[idGas],
126 ),
127 valid
128 );
129
130 label i = 0;
132 {
133 p.YGas_[j] = YGas[i]/(max(p.Y()[GAS], SMALL));
134 ++i;
135 }
136 }
137 // Populate YLiquid for each parcel
138 forAll(liquidNames, j)
139 {
140 IOField<scalar> YLiquid
141 (
142 c.fieldIOobject
143 (
144 "Y" + liquidNames[j] + stateLabels[idLiquid],
146 ),
147 valid
148 );
149
150 label i = 0;
152 {
153 p.YLiquid_[j] = YLiquid[i]/(max(p.Y()[LIQ], SMALL));
154 ++i;
155 }
156 }
157 // Populate YSolid for each parcel
159 {
160 IOField<scalar> YSolid
161 (
162 c.fieldIOobject
163 (
164 "Y" + solidNames[j] + stateLabels[idSolid],
166 ),
167 valid
168 );
169
170 label i = 0;
172 {
173 p.YSolid_[j] = YSolid[i]/(max(p.Y()[SLD], SMALL));
174 ++i;
175 }
176 }
177}
178
179
180template<class ParcelType>
181template<class CloudType>
183{
185}
186
187
188template<class ParcelType>
189template<class CloudType, class CompositionType>
191(
192 const CloudType& c,
193 const CompositionType& compModel
194)
195{
196 ParcelType::writeFields(c, compModel);
197
198 label np = c.size();
199
200 // Write the composition fractions
201 {
202 const wordList& stateLabels = compModel.stateLabels();
203
204 const label idGas = compModel.idGas();
205 const wordList& gasNames = compModel.componentNames(idGas);
206 forAll(gasNames, j)
207 {
208 IOField<scalar> YGas
209 (
210 c.fieldIOobject
211 (
212 "Y" + gasNames[j] + stateLabels[idGas],
214 ),
215 np
216 );
217
218 label i = 0;
220 {
221 YGas[i] = p0.YGas()[j]*max(p0.Y()[GAS], SMALL);
222 ++i;
223 }
224
225 YGas.write(np > 0);
226 }
227
228 const label idLiquid = compModel.idLiquid();
229 const wordList& liquidNames = compModel.componentNames(idLiquid);
230 forAll(liquidNames, j)
231 {
232 IOField<scalar> YLiquid
233 (
234 c.fieldIOobject
235 (
236 "Y" + liquidNames[j] + stateLabels[idLiquid],
238 ),
239 np
240 );
241
242 label i = 0;
244 {
245 YLiquid[i] = p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL);
246 ++i;
247 }
248
249 YLiquid.write(np > 0);
250 }
251
252 const label idSolid = compModel.idSolid();
253 const wordList& solidNames = compModel.componentNames(idSolid);
255 {
256 IOField<scalar> YSolid
257 (
258 c.fieldIOobject
259 (
260 "Y" + solidNames[j] + stateLabels[idSolid],
262 ),
263 np
264 );
265
266 label i = 0;
268 {
269 YSolid[i] = p0.YSolid()[j]*max(p0.Y()[SLD], SMALL);
270 ++i;
271 }
272
273 YSolid.write(np > 0);
274 }
275 }
276}
277
278
279template<class ParcelType>
281(
282 Ostream& os,
283 const wordRes& filters,
284 const word& delim,
285 const bool namesOnly
286) const
287{
288 ParcelType::writeProperties(os, filters, delim, namesOnly);
289
290 #undef writeProp
291 #define writeProp(Name, Value) \
292 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
293
294 writeProp("YGas", YGas_);
295 writeProp("YLiquid", YLiquid_);
296 writeProp("YSolid", YSolid_);
297 writeProp("canCombust", canCombust_);
298
299 #undef writeProp
300}
301
302
303template<class ParcelType>
304template<class CloudType>
306(
307 CloudType& c,
308 const objectRegistry& obr
309)
310{
312}
313
314
315template<class ParcelType>
316template<class CloudType>
318(
319 const CloudType& c,
320 objectRegistry& obr
321)
322{
324}
325
326
327template<class ParcelType>
328template<class CloudType, class CompositionType>
330(
331 CloudType& c,
332 const CompositionType& compModel,
333 const objectRegistry& obr
334)
335{
337
338 const label np = c.size();
339
340 // The composition fractions
341 if (np > 0)
342 {
343 const wordList& stateLabels = compModel.stateLabels();
344
345 const label idGas = compModel.idGas();
346 const wordList& gasNames = compModel.componentNames(idGas);
347 forAll(gasNames, j)
348 {
349 const word fieldName = "Y" + gasNames[j] + stateLabels[idGas];
350 const auto& YGas = cloud::lookupIOField<scalar>(fieldName, obr);
351
352 label i = 0;
354 {
355 p0.YGas()[j]*max(p0.Y()[GAS], SMALL) = YGas[i];
356 ++i;
357 }
358 }
359
360 const label idLiquid = compModel.idLiquid();
361 const wordList& liquidNames = compModel.componentNames(idLiquid);
362 forAll(liquidNames, j)
363 {
364 const word fieldName = "Y" + liquidNames[j] + stateLabels[idLiquid];
365 const auto& YLiquid = cloud::lookupIOField<scalar>(fieldName, obr);
366
367 label i = 0;
369 {
370 p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL) = YLiquid[i];
371 ++i;
372 }
373 }
374
375 const label idSolid = compModel.idSolid();
376 const wordList& solidNames = compModel.componentNames(idSolid);
378 {
379 const word fieldName = "Y" + solidNames[j] + stateLabels[idSolid];
380 const auto& YSolid = cloud::lookupIOField<scalar>(fieldName, obr);
381
382 label i = 0;
384 {
385 p0.YSolid()[j]*max(p0.Y()[SLD], SMALL) = YSolid[i];
386 ++i;
387 }
388 }
389 }
390}
391
392
393template<class ParcelType>
394template<class CloudType, class CompositionType>
396(
397 const CloudType& c,
398 const CompositionType& compModel,
399 objectRegistry& obr
400)
401{
403
404 const label np = c.size();
405
406 // Write the composition fractions
407 if (np > 0)
408 {
409 const wordList& stateLabels = compModel.stateLabels();
410
411 const label idGas = compModel.idGas();
412 const wordList& gasNames = compModel.componentNames(idGas);
413 forAll(gasNames, j)
414 {
415 const word fieldName = "Y" + gasNames[j] + stateLabels[idGas];
416 auto& YGas = cloud::createIOField<scalar>(fieldName, np, obr);
417
418 label i = 0;
420 {
421 YGas[i] = p0.YGas()[j]*max(p0.Y()[GAS], SMALL);
422 ++i;
423 }
424 }
425
426 const label idLiquid = compModel.idLiquid();
427 const wordList& liquidNames = compModel.componentNames(idLiquid);
428 forAll(liquidNames, j)
429 {
430 const word fieldName = "Y" + liquidNames[j] + stateLabels[idLiquid];
431 auto& YLiquid = cloud::createIOField<scalar>(fieldName, np, obr);
432
433 label i = 0;
435 {
436 YLiquid[i] = p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL);
437 ++i;
438 }
439 }
440
441 const label idSolid = compModel.idSolid();
442 const wordList& solidNames = compModel.componentNames(idSolid);
444 {
445 const word fieldName = "Y" + solidNames[j] + stateLabels[idSolid];
446 auto& YSolid = cloud::createIOField<scalar>(fieldName, np, obr);
447
448 label i = 0;
450 {
451 YSolid[i] = p0.YSolid()[j]*max(p0.Y()[SLD], SMALL);
452 ++i;
453 }
454 }
455 }
456}
457
458
459// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
460
461template<class ParcelType>
462Foam::Ostream& Foam::operator<<
463(
464 Ostream& os,
466)
467{
468 scalarField YGasLoc(p.YGas());
469 scalarField YLiquidLoc(p.YLiquid());
470 scalarField YSolidLoc(p.YSolid());
471
472 if (os.format() == IOstream::ASCII)
473 {
474 os << static_cast<const ParcelType&>(p)
475 << token::SPACE << YGasLoc
476 << token::SPACE << YLiquidLoc
477 << token::SPACE << YSolidLoc;
478 }
479 else
480 {
481 os << static_cast<const ParcelType&>(p);
482 os << YGasLoc << YLiquidLoc << YSolidLoc;
483 }
484
486 return os;
487}
488
489
490// ************************************************************************* //
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.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase.
scalarField YLiquid_
Mass fractions of liquids [].
scalarField YSolid_
Mass fractions of solids [].
static const std::size_t sizeofFields
Size in bytes of the fields.
scalarField YGas_
Mass fractions of gases [].
static void readObjects(CloudType &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=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
const volScalarField & p0
Definition: EEqn.H:36
dynamicFvMesh & mesh
const wordList solidNames(rp["solid"])
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define FUNCTION_NAME
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
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