SprayParcelIO.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "SprayParcel.H"
27 #include "IOstreams.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<class ParcelType>
33 (
34  sizeof(SprayParcel<ParcelType>) - sizeof(ParcelType)
35 );
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class ParcelType>
42 (
43  const polyMesh& mesh,
44  Istream& is,
45  bool readFields
46 )
47 :
48  ParcelType(mesh, is, readFields),
49  d0_(0.0),
50  position0_(Zero),
51  sigma_(0.0),
52  mu_(0.0),
53  liquidCore_(0.0),
54  KHindex_(0.0),
55  y_(0.0),
56  yDot_(0.0),
57  tc_(0.0),
58  ms_(0.0),
59  injector_(1.0),
60  tMom_(great),
61  user_(0.0)
62 {
63  if (readFields)
64  {
65  if (is.format() == IOstream::ASCII)
66  {
67  d0_ = readScalar(is);
68  is >> position0_;
69  sigma_ = readScalar(is);
70  mu_ = readScalar(is);
71  liquidCore_ = readScalar(is);
72  KHindex_ = readScalar(is);
73  y_ = readScalar(is);
74  yDot_ = readScalar(is);
75  tc_ = readScalar(is);
76  ms_ = readScalar(is);
77  injector_ = readScalar(is);
78  tMom_ = readScalar(is);
79  user_ = readScalar(is);
80  }
81  else
82  {
83  is.read(reinterpret_cast<char*>(&d0_), sizeofFields_);
84  }
85  }
86 
87  // Check state of Istream
88  is.check
89  (
90  "SprayParcel<ParcelType>::SprayParcel"
91  "("
92  "const polyMesh, "
93  "Istream&, "
94  "bool"
95  ")"
96  );
97 }
98 
99 
100 template<class ParcelType>
101 template<class CloudType>
103 {
105 }
106 
107 
108 template<class ParcelType>
109 template<class CloudType, class CompositionType>
111 (
112  CloudType& c,
113  const CompositionType& compModel
114 )
115 {
116  bool write = c.size();
117 
118  ParcelType::readFields(c, compModel);
119 
120  IOField<scalar> d0(c.fieldIOobject("d0", IOobject::MUST_READ), write);
121  c.checkFieldIOobject(c, d0);
122 
123  IOField<vector> position0
124  (
125  c.fieldIOobject("position0", IOobject::MUST_READ),
126  write
127  );
128  c.checkFieldIOobject(c, position0);
129 
130  IOField<scalar> sigma(c.fieldIOobject("sigma", IOobject::MUST_READ), write);
132 
133  IOField<scalar> mu(c.fieldIOobject("mu", IOobject::MUST_READ), write);
134  c.checkFieldIOobject(c, mu);
135 
136  IOField<scalar> liquidCore
137  (
138  c.fieldIOobject("liquidCore", IOobject::MUST_READ),
139  write
140  );
141  c.checkFieldIOobject(c, liquidCore);
142 
143  IOField<scalar> KHindex
144  (
145  c.fieldIOobject("KHindex", IOobject::MUST_READ),
146  write
147  );
148  c.checkFieldIOobject(c, KHindex);
149 
151  (
152  c.fieldIOobject("y", IOobject::MUST_READ),
153  write
154  );
155  c.checkFieldIOobject(c, y);
156 
157  IOField<scalar> yDot
158  (
159  c.fieldIOobject("yDot", IOobject::MUST_READ),
160  write
161  );
162  c.checkFieldIOobject(c, yDot);
163 
164  IOField<scalar> tc
165  (
166  c.fieldIOobject("tc", IOobject::MUST_READ),
167  write
168  );
169  c.checkFieldIOobject(c, tc);
170 
171  IOField<scalar> ms
172  (
173  c.fieldIOobject("ms", IOobject::MUST_READ),
174  write
175  );
176  c.checkFieldIOobject(c, ms);
177 
178  IOField<scalar> injector
179  (
180  c.fieldIOobject("injector", IOobject::MUST_READ),
181  write
182  );
183  c.checkFieldIOobject(c, injector);
184 
185  IOField<scalar> tMom
186  (
187  c.fieldIOobject("tMom", IOobject::MUST_READ),
188  write
189  );
190  c.checkFieldIOobject(c, tMom);
191 
192  IOField<scalar> user
193  (
194  c.fieldIOobject("user", IOobject::MUST_READ),
195  write
196  );
197  c.checkFieldIOobject(c, user);
198 
199  label i = 0;
200  forAllIter(typename CloudType, c, iter)
201  {
202  SprayParcel<ParcelType>& p = iter();
203  p.d0_ = d0[i];
204  p.position0_ = position0[i];
205  p.sigma_ = sigma[i];
206  p.mu_ = mu[i];
207  p.liquidCore_ = liquidCore[i];
208  p.KHindex_ = KHindex[i];
209  p.y_ = y[i];
210  p.yDot_ = yDot[i];
211  p.tc_ = tc[i];
212  p.ms_ = ms[i];
213  p.injector_ = injector[i];
214  p.tMom_ = tMom[i];
215  p.user_ = user[i];
216  i++;
217  }
218 }
219 
220 
221 template<class ParcelType>
222 template<class CloudType>
224 {
225  ParcelType::writeFields(c);
226 }
227 
228 
229 template<class ParcelType>
230 template<class CloudType, class CompositionType>
232 (
233  const CloudType& c,
234  const CompositionType& compModel
235 )
236 {
237  ParcelType::writeFields(c, compModel);
238 
239  label np = c.size();
240 
241  IOField<scalar> d0(c.fieldIOobject("d0", IOobject::NO_READ), np);
242  IOField<vector> position0
243  (
244  c.fieldIOobject("position0", IOobject::NO_READ),
245  np
246  );
247  IOField<scalar> sigma(c.fieldIOobject("sigma", IOobject::NO_READ), np);
248  IOField<scalar> mu(c.fieldIOobject("mu", IOobject::NO_READ), np);
249  IOField<scalar> liquidCore
250  (
251  c.fieldIOobject("liquidCore", IOobject::NO_READ),
252  np
253  );
254  IOField<scalar> KHindex(c.fieldIOobject("KHindex", IOobject::NO_READ), np);
255  IOField<scalar> y(c.fieldIOobject("y", IOobject::NO_READ), np);
256  IOField<scalar> yDot(c.fieldIOobject("yDot", IOobject::NO_READ), np);
257  IOField<scalar> tc(c.fieldIOobject("tc", IOobject::NO_READ), np);
258  IOField<scalar> ms(c.fieldIOobject("ms", IOobject::NO_READ), np);
259  IOField<scalar> injector
260  (
261  c.fieldIOobject("injector", IOobject::NO_READ),
262  np
263  );
264  IOField<scalar> tMom(c.fieldIOobject("tMom", IOobject::NO_READ), np);
265  IOField<scalar> user(c.fieldIOobject("user", IOobject::NO_READ), np);
266 
267  label i = 0;
268  forAllConstIter(typename CloudType, c, iter)
269  {
270  const SprayParcel<ParcelType>& p = iter();
271  d0[i] = p.d0_;
272  position0[i] = p.position0_;
273  sigma[i] = p.sigma_;
274  mu[i] = p.mu_;
275  liquidCore[i] = p.liquidCore_;
276  KHindex[i] = p.KHindex_;
277  y[i] = p.y_;
278  yDot[i] = p.yDot_;
279  tc[i] = p.tc_;
280  ms[i] = p.ms_;
281  injector[i] = p.injector_;
282  tMom[i] = p.tMom_;
283  user[i] = p.user_;
284  i++;
285  }
286 
287  const bool write = np > 0;
288 
289  d0.write(write);
290  position0.write(write);
291  sigma.write(write);
292  mu.write(write);
293  liquidCore.write(write);
294  KHindex.write(write);
295  y.write(write);
296  yDot.write(write);
297  tc.write(write);
298  ms.write(write);
299  injector.write(write);
300  tMom.write(write);
301  user.write(write);
302 }
303 
304 
305 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
306 
307 template<class ParcelType>
308 Foam::Ostream& Foam::operator<<
309 (
310  Ostream& os,
312 )
313 {
314  if (os.format() == IOstream::ASCII)
315  {
316  os << static_cast<const ParcelType&>(p)
317  << token::SPACE << p.d0()
318  << token::SPACE << p.position0()
319  << token::SPACE << p.sigma()
320  << token::SPACE << p.mu()
321  << token::SPACE << p.liquidCore()
322  << token::SPACE << p.KHindex()
323  << token::SPACE << p.y()
324  << token::SPACE << p.yDot()
325  << token::SPACE << p.tc()
326  << token::SPACE << p.ms()
327  << token::SPACE << p.injector()
328  << token::SPACE << p.tMom()
329  << token::SPACE << p.user();
330  }
331  else
332  {
333  os << static_cast<const ParcelType&>(p);
334  os.write
335  (
336  reinterpret_cast<const char*>(&p.d0_),
338  );
339  }
340 
341  // Check state of Ostream
342  os.check
343  (
344  "Ostream& operator<<(Ostream&, const SprayParcel<ParcelType>&)"
345  );
346 
347  return os;
348 }
349 
350 
351 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.)
Definition: SprayParcel.H:183
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:148
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet )
Definition: SprayParcel.H:161
SprayParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: SprayParcelI.H:108
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
scalar y
virtual Istream & read(token &)=0
Return next token from stream.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
static const zero Zero
Definition: zero.H:97
scalar injector_
Injected from injector (needed e.g. for calculating distance.
Definition: SprayParcel.H:180
streamFormat format() const
Return current stream format.
Definition: IOstream.H:374
scalar user_
Passive scalar (extra variable to be defined by user)
Definition: SprayParcel.H:186
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const dimensionedScalar mu
Atomic mass unit.
scalar sigma_
Liquid surface tension [N/m].
Definition: SprayParcel.H:155
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
vector position0_
Injection position.
Definition: SprayParcel.H:152
scalar yDot_
Rate of change of spherical deviation.
Definition: SprayParcel.H:170
scalar ms_
Stripped parcel mass due to breakup.
Definition: SprayParcel.H:176
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:187
scalar d0_
Initial droplet diameter.
Definition: SprayParcel.H:149
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar tc_
Characteristic time (used in atomisation and/or breakup model)
Definition: SprayParcel.H:173
volScalarField & p
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:158
scalar y_
Spherical deviation.
Definition: SprayParcel.H:167
scalar KHindex_
Index for KH Breakup.
Definition: SprayParcel.H:164
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:167
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Reaching spray parcel, with added functionality for atomisation and breakup.
Definition: SprayParcel.H:43