Loading...
Searching...
No Matches
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-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
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),
62{
63 if (readFields)
64 {
68
69 is >> Yg >> Yl >> Ys;
70
71 YGas_.transfer(Yg);
72 YLiquid_.transfer(Yl);
73 YSolid_.transfer(Ys);
74 }
77}
78
79
80template<class ParcelType>
81template<class CloudType>
84 ParcelType::readFields(c);
85}
86
87
88template<class ParcelType>
89template<class CloudType, class CompositionType>
91(
92 CloudType& c,
93 const CompositionType& compModel
94)
95{
96 const bool readOnProc = 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.newIOobject
123 (
124 "Y" + gasNames[j] + stateLabels[idGas],
126 ),
127 readOnProc
128 );
129
130 label i = 0;
131 for (ReactingMultiphaseParcel<ParcelType>& p : c)
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.newIOobject
143 (
144 "Y" + liquidNames[j] + stateLabels[idLiquid],
146 ),
147 readOnProc
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.newIOobject
163 (
164 "Y" + solidNames[j] + stateLabels[idSolid],
166 ),
167 readOnProc
168 );
169
170 label i = 0;
172 {
173 p.YSolid_[j] = YSolid[i]/(max(p.Y()[SLD], SMALL));
174 ++i;
176 }
177}
178
179
180template<class ParcelType>
181template<class CloudType>
184 ParcelType::writeFields(c);
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 const label np = c.size();
199 const bool writeOnProc = c.size();
200
201 // Write the composition fractions
202 {
203 const wordList& stateLabels = compModel.stateLabels();
204
205 const label idGas = compModel.idGas();
206 const wordList& gasNames = compModel.componentNames(idGas);
207 forAll(gasNames, j)
208 {
209 IOField<scalar> YGas
210 (
211 c.newIOobject
212 (
213 "Y" + gasNames[j] + stateLabels[idGas],
215 ),
216 np
217 );
218
219 label i = 0;
221 {
222 YGas[i] = p0.YGas()[j]*max(p0.Y()[GAS], SMALL);
223 ++i;
224 }
225
226 YGas.write(writeOnProc);
227 }
228
229 const label idLiquid = compModel.idLiquid();
230 const wordList& liquidNames = compModel.componentNames(idLiquid);
231 forAll(liquidNames, j)
232 {
233 IOField<scalar> YLiquid
234 (
235 c.newIOobject
236 (
237 "Y" + liquidNames[j] + stateLabels[idLiquid],
239 ),
240 np
241 );
242
243 label i = 0;
245 {
246 YLiquid[i] = p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL);
247 ++i;
248 }
249
250 YLiquid.write(writeOnProc);
251 }
252
253 const label idSolid = compModel.idSolid();
254 const wordList& solidNames = compModel.componentNames(idSolid);
256 {
257 IOField<scalar> YSolid
258 (
259 c.newIOobject
260 (
261 "Y" + solidNames[j] + stateLabels[idSolid],
263 ),
264 np
265 );
266
267 label i = 0;
269 {
270 YSolid[i] = p0.YSolid()[j]*max(p0.Y()[SLD], SMALL);
271 ++i;
272 }
273
274 YSolid.write(writeOnProc);
275 }
276 }
277}
278
279
280template<class ParcelType>
282(
283 Ostream& os,
284 const wordRes& filters,
285 const word& delim,
286 const bool namesOnly
287) const
288{
289 ParcelType::writeProperties(os, filters, delim, namesOnly);
290
291 #undef writeProp
292 #define writeProp(Name, Value) \
293 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
294
295 writeProp("YGas", YGas_);
296 writeProp("YLiquid", YLiquid_);
297 writeProp("YSolid", YSolid_);
298 writeProp("canCombust", canCombust_);
300 #undef writeProp
301}
302
303
304template<class ParcelType>
305template<class CloudType>
307(
308 CloudType& c,
309 const objectRegistry& obr
310)
312 ParcelType::readObjects(c, obr);
313}
314
315
316template<class ParcelType>
317template<class CloudType>
319(
320 const CloudType& c,
321 objectRegistry& obr
322)
324 ParcelType::writeObjects(c, obr);
325}
326
327
328template<class ParcelType>
329template<class CloudType, class CompositionType>
331(
332 CloudType& c,
333 const CompositionType& compModel,
334 const objectRegistry& obr
335)
336{
337 ParcelType::readObjects(c, obr);
338
339 // const label np = c.size();
340 const bool readOnProc = c.size();
341
342 // The composition fractions
343 if (readOnProc)
344 {
345 const wordList& stateLabels = compModel.stateLabels();
346
347 const label idGas = compModel.idGas();
348 const wordList& gasNames = compModel.componentNames(idGas);
349 forAll(gasNames, j)
350 {
351 const word fieldName = "Y" + gasNames[j] + stateLabels[idGas];
352 const auto& YGas = cloud::lookupIOField<scalar>(fieldName, obr);
353
354 label i = 0;
356 {
357 p0.YGas()[j]*max(p0.Y()[GAS], SMALL) = YGas[i];
358 ++i;
359 }
360 }
361
362 const label idLiquid = compModel.idLiquid();
363 const wordList& liquidNames = compModel.componentNames(idLiquid);
364 forAll(liquidNames, j)
365 {
366 const word fieldName = "Y" + liquidNames[j] + stateLabels[idLiquid];
367 const auto& YLiquid = cloud::lookupIOField<scalar>(fieldName, obr);
368
369 label i = 0;
371 {
372 p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL) = YLiquid[i];
373 ++i;
374 }
375 }
376
377 const label idSolid = compModel.idSolid();
378 const wordList& solidNames = compModel.componentNames(idSolid);
380 {
381 const word fieldName = "Y" + solidNames[j] + stateLabels[idSolid];
382 const auto& YSolid = cloud::lookupIOField<scalar>(fieldName, obr);
383
384 label i = 0;
386 {
387 p0.YSolid()[j]*max(p0.Y()[SLD], SMALL) = YSolid[i];
388 ++i;
389 }
391 }
392}
393
394
395template<class ParcelType>
396template<class CloudType, class CompositionType>
398(
399 const CloudType& c,
400 const CompositionType& compModel,
401 objectRegistry& obr
402)
403{
404 ParcelType::writeObjects(c, obr);
405
406 const label np = c.size();
407 const bool writeOnProc = c.size();
408
409 // Write the composition fractions
410 if (writeOnProc)
411 {
412 const wordList& stateLabels = compModel.stateLabels();
413
414 const label idGas = compModel.idGas();
415 const wordList& gasNames = compModel.componentNames(idGas);
416 forAll(gasNames, j)
417 {
418 const word fieldName = "Y" + gasNames[j] + stateLabels[idGas];
419 auto& YGas = cloud::createIOField<scalar>(fieldName, np, obr);
420
421 label i = 0;
423 {
424 YGas[i] = p0.YGas()[j]*max(p0.Y()[GAS], SMALL);
425 ++i;
426 }
427 }
428
429 const label idLiquid = compModel.idLiquid();
430 const wordList& liquidNames = compModel.componentNames(idLiquid);
431 forAll(liquidNames, j)
432 {
433 const word fieldName = "Y" + liquidNames[j] + stateLabels[idLiquid];
434 auto& YLiquid = cloud::createIOField<scalar>(fieldName, np, obr);
435
436 label i = 0;
438 {
439 YLiquid[i] = p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL);
440 ++i;
441 }
442 }
443
444 const label idSolid = compModel.idSolid();
445 const wordList& solidNames = compModel.componentNames(idSolid);
447 {
448 const word fieldName = "Y" + solidNames[j] + stateLabels[idSolid];
449 auto& YSolid = cloud::createIOField<scalar>(fieldName, np, obr);
450
451 label i = 0;
453 {
454 YSolid[i] = p0.YSolid()[j]*max(p0.Y()[SLD], SMALL);
455 ++i;
456 }
457 }
458 }
459}
460
461
462// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
463
464template<class ParcelType>
465Foam::Ostream& Foam::operator<<
466(
467 Ostream& os,
469)
470{
471 scalarField YGasLoc(p.YGas());
472 scalarField YLiquidLoc(p.YLiquid());
473 scalarField YSolidLoc(p.YSolid());
474
476 {
478 << token::SPACE << YGasLoc
479 << token::SPACE << YLiquidLoc
480 << token::SPACE << YSolidLoc;
481 }
482 else
483 {
485 os << YGasLoc << YLiquidLoc << YSolidLoc;
486 }
487
489 return os;
490}
491
492
493// ************************************************************************* //
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.
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
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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 [].
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.
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 YGas_
Mass fractions of gases [].
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
const scalarField & YGas() const
Return const access to mass fractions of gases.
ReactingMultiphaseParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
const scalarField & YSolid() const
Return const access to mass fractions of solids.
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
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
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
const dimensionedScalar c
Speed of light in a vacuum.
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DSMCCloud< dsmcParcel > CloudType
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299