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 -------------------------------------------------------------------------------
10 License
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 
34 template<class ParcelType>
37 
38 
39 template<class ParcelType>
41 (
42  sizeof(scalar)
43 );
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 template<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),
59  canCombust_(1)
60 {
61  if (readFields)
62  {
64 
65  is >> F;
66 
67  F_.transfer(F);
68  }
69 
70  is.check(FUNCTION_NAME);
71 }
72 
73 
74 template<class ParcelType>
75 template<class CloudType>
77 {
79 }
80 
81 
82 template<class ParcelType>
83 template<class CloudType, class CompositionType>
85 (
86  CloudType& c,
87  const CompositionType& compModel
88 )
89 {
90  const bool readOnProc = c.size();
91 
93 
94  IOField<scalar> mass0
95  (
96  c.fieldIOobject("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
124  IOField<scalar> F
125  (
126  c.fieldIOobject
127  (
128  "F" + name(i),
130  ),
131  readOnProc
132  );
133 
134  label j = 0;
135  for (ReactingHeterogeneousParcel<ParcelType>& p : c)
136  {
137  p.F_[i] = F[j];
138  ++j;
139  }
140  }
141 
142 
143  forAll(solidNames, j)
144  {
145  IOField<scalar> Y
146  (
147  c.fieldIOobject
148  (
149  "Y" + solidNames[j],
151  ),
152  readOnProc
153  );
154 
155  label i = 0;
156  for (ReactingHeterogeneousParcel<ParcelType>& p : c)
157  {
158  p.Y()[j] = Y[i];
159  ++i;
160  }
161  }
162 }
163 
164 
165 template<class ParcelType>
166 template<class CloudType>
168 (
169  const CloudType& c
170 )
171 {
173 }
174 
175 
176 template<class ParcelType>
177 template<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.fieldIOobject("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  {
209  IOField<scalar> F
210  (
211  c.fieldIOobject
212  (
213  "F" + name(i),
215  ),
216  np
217  );
218 
219  for (const ReactingHeterogeneousParcel<ParcelType>& p0 : c)
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 
230  forAll(solidNames, j)
231  {
232  IOField<scalar> Y
233  (
234  c.fieldIOobject
235  (
236  "Y" + solidNames[j],
238  ),
239  np
240  );
241 
242  label i = 0;
243  for (const ReactingHeterogeneousParcel<ParcelType>& p0 : c)
244  {
245  Y[i] = p0.Y()[j];
246  ++i;
247  }
248 
249  Y.write(writeOnProc);
250  }
251 }
252 
253 
254 template<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_);
271 
272  #undef writeProp
273 }
274 
275 
276 template<class ParcelType>
277 template<class CloudType>
279 (
280  CloudType& c,
281  const objectRegistry& obr
282 )
283 {
284  ParcelType::readObjects(c, obr);
285 }
286 
287 
288 template<class ParcelType>
289 template<class CloudType>
291 (
292  const CloudType& c,
293  objectRegistry& obr
294 )
295 {
296  ParcelType::writeObjects(c, obr);
297 }
298 
299 
300 template<class ParcelType>
301 template<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;
318 }
319 
320 
321 template<class ParcelType>
322 template<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);
362  forAll(solidNames, j)
363  {
364  const word fieldName = "Y" + solidNames[j];
365  auto& Y = cloud::createIOField<scalar>(fieldName, np, obr);
366 
367  label i = 0;
368  for (const ReactingHeterogeneousParcel<ParcelType>& p0 : c)
369  {
370  Y[i] = p0.Y()[j];
371  ++i;
372  }
373  }
374  }
375 }
376 
377 
378 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
379 
380 template<class ParcelType>
381 Foam::Ostream& Foam::operator<<
382 (
383  Ostream& os,
384  const ReactingHeterogeneousParcel<ParcelType>& p
385 )
386 {
387  scalarField F(p.F());
389  {
390  os << static_cast<const ParcelType&>(p)
391  << token::SPACE << F;
392  }
393  else
394  {
395  os << static_cast<const ParcelType&>(p);
396  os << F ;
397  }
398 
400  return os;
401 }
402 
403 
404 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
DSMCCloud< dsmcParcel > CloudType
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
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:57
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
"ascii" (normal default)
static const std::size_t sizeofFields
Size in bytes of the fields.
#define writeProp(Name, Value)
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
const wordList solidNames(rp["solid"])
ReactingHeterogeneousParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
const scalarField & F() const
Return const access to F.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Space [isspace].
Definition: token.H:131
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
volVectorField F(fluid.F())
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
List< word > wordList
List of word.
Definition: fileName.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
PtrList< volScalarField > & Y
const dimensionedScalar c
Speed of light in a vacuum.
Nothing to be read.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
scalarField F_
Progress variables for reactions.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
Registry of regIOobjects.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
A class for handling character strings derived from std::string.
Definition: string.H:72
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
A primitive field of type <T> with automated input and output.
streamFormat format() const noexcept
Get the current stream format.
const volScalarField & p0
Definition: EEqn.H:36
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
Definition: ThermoParcel.H:53
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127