Loading...
Searching...
No Matches
ReactingHeterogeneousParcelIO.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) 2018-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#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 F_(0),
60{
61 if (readFields)
62 {
64
65 is >> F;
66
67 F_.transfer(F);
68 }
71}
72
73
74template<class ParcelType>
75template<class CloudType>
78 ParcelType::readFields(c);
79}
80
81
82template<class ParcelType>
83template<class CloudType, class CompositionType>
85(
86 CloudType& c,
87 const CompositionType& compModel
88)
89{
90 const bool readOnProc = c.size();
91
92 ParcelType::readFields(c);
93
94 IOField<scalar> mass0
95 (
96 c.newIOobject("mass0", IOobject::MUST_READ),
97 readOnProc
98 );
99 c.checkFieldIOobject(c, mass0);
100
101 label i = 0;
103 {
104 p.mass0() = mass0[i];
105 ++i;
106 }
107
108 const label idSolid = compModel.idSolid();
109 const wordList& solidNames = compModel.componentNames(idSolid);
110
111 // WIP until find a way to get info from Reacting Heterogeneous model
112 label nF(1);
113
114 // Set storage for each Y... for each parcel
115 for (ReactingHeterogeneousParcel<ParcelType>& p : c)
116 {
117 p.Y().setSize(solidNames.size(), Zero);
118 p.F_.setSize(nF, Zero);
119 }
120
121 for (label i = 0; i < nF; i++)
122 {
123 // Read F
125 (
126 c.newIOobject
127 (
128 "F" + name(i),
130 ),
131 readOnProc
132 );
133
134 label j = 0;
136 {
137 p.F_[i] = F[j];
138 ++j;
139 }
140 }
141
142
144 {
146 (
147 c.newIOobject
148 (
149 "Y" + solidNames[j],
151 ),
152 readOnProc
153 );
154
155 label i = 0;
157 {
158 p.Y()[j] = Y[i];
159 ++i;
161 }
162}
163
164
165template<class ParcelType>
166template<class CloudType>
168(
169 const CloudType& c
170)
172 ParcelType::writeFields(c);
173}
174
175
176template<class ParcelType>
177template<class CloudType, class CompositionType>
179(
180 const CloudType& c,
181 const CompositionType& compModel
182)
183{
184 // Skip Reacting layer. This class takes charge of
185 // writing Ysolids and F
187
188 const label np = c.size();
189 const bool writeOnProc = c.size();
190
191 IOField<scalar> mass0(c.newIOobject("mass0", IOobject::NO_READ), np);
192
193 label nF = 0;
194 label i = 0;
196 {
197 mass0[i] = p.mass0_;
198 if (!i)
199 {
200 // Assume same size throughout
201 nF = p.F().size();
202 }
203 ++i;
204 }
205 mass0.write(writeOnProc);
206
207 for (label i = 0; i < nF; i++)
208 {
210 (
211 c.newIOobject
212 (
213 "F" + name(i),
215 ),
216 np
217 );
218
220 {
221 F = p0.F()[i];
222 }
223
224 F.write(writeOnProc);
225 }
226
227 const label idSolid = compModel.idSolid();
228 const wordList& solidNames = compModel.componentNames(idSolid);
229
231 {
233 (
234 c.newIOobject
235 (
236 "Y" + solidNames[j],
238 ),
239 np
240 );
241
242 label i = 0;
244 {
245 Y[i] = p0.Y()[j];
246 ++i;
247 }
249 Y.write(writeOnProc);
250 }
251}
252
253
254template<class ParcelType>
256(
257 Ostream& os,
258 const wordRes& filters,
259 const word& delim,
260 const bool namesOnly
261) const
262{
263 ParcelType::writeProperties(os, filters, delim, namesOnly);
264
265 #undef writeProp
266 #define writeProp(Name, Value) \
267 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
268
269 writeProp("F", F_);
270 writeProp("canCombust", canCombust_);
272 #undef writeProp
273}
274
275
276template<class ParcelType>
277template<class CloudType>
279(
280 CloudType& c,
281 const objectRegistry& obr
282)
284 ParcelType::readObjects(c, obr);
285}
286
287
288template<class ParcelType>
289template<class CloudType>
291(
292 const CloudType& c,
293 objectRegistry& obr
294)
296 ParcelType::writeObjects(c, obr);
297}
298
299
300template<class ParcelType>
301template<class CloudType, class CompositionType>
303(
304 CloudType& c,
305 const CompositionType& compModel,
306 const objectRegistry& obr
307)
308{
309 //ParcelType::readObjects(c, obr);
310 // Skip Reacting layer
311 ThermoParcel<KinematicParcel<particle>>::readObjects(c, obr);
312
313 // const bool readOnProc = c.size();
314
316 << "Reading of objects is still a work-in-progress" << nl;
317
318}
319
320
321template<class ParcelType>
322template<class CloudType, class CompositionType>
324(
325 const CloudType& c,
326 const CompositionType& compModel,
327 objectRegistry& obr
328)
329{
330 //ParcelType::writeObjects(c, obr);
331 // Skip Reacting layer
332 ThermoParcel<KinematicParcel<particle>>::writeObjects(c, obr);
333
334 const label np = c.size();
335 const bool writeOnProc = c.size();
336
337 // WIP
338 label nF = 0;
340 {
341 nF = p0.F().size();
342 break;
343 }
344
345 if (writeOnProc)
346 {
347 for (label i = 0; i < nF; i++)
348 {
349 const word fieldName = "F" + name(i);
350 auto& F = cloud::createIOField<scalar>(fieldName, np, obr);
351
352 label j = 0;
353 for (const ReactingHeterogeneousParcel<ParcelType>& p0 : c)
354 {
355 F[j] = p0.F()[i];
356 ++j;
357 }
358 }
359
360 const label idSolid = compModel.idSolid();
361 const wordList& solidNames = compModel.componentNames(idSolid);
363 {
364 const word fieldName = "Y" + solidNames[j];
365 auto& Y = cloud::createIOField<scalar>(fieldName, np, obr);
366
367 label i = 0;
369 {
370 Y[i] = p0.Y()[j];
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 scalarField F(p.F());
389 {
391 << token::SPACE << F;
392 }
393 else
394 {
396 os << F ;
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...
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
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.
Definition polyMesh.H:79
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
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
const volScalarField & p0
Definition EEqn.H:36
dynamicFvMesh & mesh
const wordList solidNames(rp["solid"])
OBJstream os(runTime.globalPath()/outputName)
auto & name
#define writeProp(Name, Value)
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
volVectorField F(fluid.F())
const dimensionedScalar c
Speed of light in a vacuum.
List< word > wordList
List of word.
Definition fileName.H:60
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).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299