SprayParcelIO.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 | Copyright (C) 2011-2015 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>
34 
35 
36 template<class ParcelType>
38 (
39  sizeof(SprayParcel<ParcelType>) - sizeof(ParcelType)
40 );
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 template<class ParcelType>
47 (
48  const polyMesh& mesh,
49  Istream& is,
50  bool readFields
51 )
52 :
53  ParcelType(mesh, is, readFields),
54  d0_(0.0),
55  position0_(vector::zero),
56  sigma_(0.0),
57  mu_(0.0),
58  liquidCore_(0.0),
59  KHindex_(0.0),
60  y_(0.0),
61  yDot_(0.0),
62  tc_(0.0),
63  ms_(0.0),
64  injector_(1.0),
65  tMom_(GREAT),
66  user_(0.0)
67 {
68  if (readFields)
69  {
70  if (is.format() == IOstream::ASCII)
71  {
72  d0_ = readScalar(is);
73  is >> position0_;
74  sigma_ = readScalar(is);
75  mu_ = readScalar(is);
76  liquidCore_ = readScalar(is);
77  KHindex_ = readScalar(is);
78  y_ = readScalar(is);
79  yDot_ = readScalar(is);
80  tc_ = readScalar(is);
81  ms_ = readScalar(is);
82  injector_ = readScalar(is);
83  tMom_ = readScalar(is);
84  user_ = readScalar(is);
85  }
86  else
87  {
88  is.read(reinterpret_cast<char*>(&d0_), sizeofFields_);
89  }
90  }
91 
92  // Check state of Istream
93  is.check
94  (
95  "SprayParcel<ParcelType>::SprayParcel"
96  "("
97  "const polyMesh&, "
98  "Istream&, "
99  "bool"
100  ")"
101  );
102 }
103 
104 
105 template<class ParcelType>
106 template<class CloudType>
108 {
109  if (!c.size())
110  {
111  return;
112  }
113 
115 }
116 
117 
118 template<class ParcelType>
119 template<class CloudType, class CompositionType>
121 (
122  CloudType& c,
123  const CompositionType& compModel
124 )
125 {
126  if (!c.size())
127  {
128  return;
129  }
130 
131  ParcelType::readFields(c, compModel);
132 
133  IOField<scalar> d0(c.fieldIOobject("d0", IOobject::MUST_READ));
134  c.checkFieldIOobject(c, d0);
135 
136  IOField<vector> position0
137  (
138  c.fieldIOobject("position0", IOobject::MUST_READ)
139  );
140  c.checkFieldIOobject(c, position0);
141 
142  IOField<scalar> sigma(c.fieldIOobject("sigma", IOobject::MUST_READ));
143  c.checkFieldIOobject(c, sigma);
144 
145  IOField<scalar> mu(c.fieldIOobject("mu", IOobject::MUST_READ));
146  c.checkFieldIOobject(c, mu);
147 
148  IOField<scalar> liquidCore(c.fieldIOobject
149  (
150  "liquidCore", IOobject::MUST_READ)
151  );
152  c.checkFieldIOobject(c, liquidCore);
153 
154  IOField<scalar> KHindex(c.fieldIOobject("KHindex", IOobject::MUST_READ));
155  c.checkFieldIOobject(c, KHindex);
156 
157  IOField<scalar> y(c.fieldIOobject("y", IOobject::MUST_READ));
158  c.checkFieldIOobject(c, y);
159 
160  IOField<scalar> yDot(c.fieldIOobject("yDot", IOobject::MUST_READ));
161  c.checkFieldIOobject(c, yDot);
162 
163  IOField<scalar> tc(c.fieldIOobject("tc", IOobject::MUST_READ));
164  c.checkFieldIOobject(c, tc);
165 
166  IOField<scalar> ms(c.fieldIOobject("ms", IOobject::MUST_READ));
167  c.checkFieldIOobject(c, ms);
168 
169  IOField<scalar> injector(c.fieldIOobject("injector", IOobject::MUST_READ));
170  c.checkFieldIOobject(c, injector);
171 
172  IOField<scalar> tMom(c.fieldIOobject("tMom", IOobject::MUST_READ));
173  c.checkFieldIOobject(c, tMom);
174 
175  IOField<scalar> user(c.fieldIOobject("user", IOobject::MUST_READ));
176  c.checkFieldIOobject(c, user);
177 
178  label i = 0;
179  forAllIter(typename Cloud<SprayParcel<ParcelType> >, c, iter)
180  {
181  SprayParcel<ParcelType>& p = iter();
182  p.d0_ = d0[i];
183  p.position0_ = position0[i];
184  p.sigma_ = sigma[i];
185  p.mu_ = mu[i];
186  p.liquidCore_ = liquidCore[i];
187  p.KHindex_ = KHindex[i];
188  p.y_ = y[i];
189  p.yDot_ = yDot[i];
190  p.tc_ = tc[i];
191  p.ms_ = ms[i];
192  p.injector_ = injector[i];
193  p.tMom_ = tMom[i];
194  p.user_ = user[i];
195  i++;
196  }
197 }
198 
199 
200 template<class ParcelType>
201 template<class CloudType>
203 {
204  ParcelType::writeFields(c);
205 }
206 
207 
208 template<class ParcelType>
209 template<class CloudType, class CompositionType>
211 (
212  const CloudType& c,
213  const CompositionType& compModel
214 )
215 {
216  ParcelType::writeFields(c, compModel);
217 
218  label np = c.size();
219 
220  IOField<scalar> d0(c.fieldIOobject("d0", IOobject::NO_READ), np);
221  IOField<vector> position0
222  (
223  c.fieldIOobject("position0", IOobject::NO_READ),
224  np
225  );
226  IOField<scalar> sigma(c.fieldIOobject("sigma", IOobject::NO_READ), np);
227  IOField<scalar> mu(c.fieldIOobject("mu", IOobject::NO_READ), np);
228  IOField<scalar> liquidCore
229  (
230  c.fieldIOobject("liquidCore", IOobject::NO_READ),
231  np
232  );
233  IOField<scalar> KHindex(c.fieldIOobject("KHindex", IOobject::NO_READ), np);
234  IOField<scalar> y(c.fieldIOobject("y", IOobject::NO_READ), np);
235  IOField<scalar> yDot(c.fieldIOobject("yDot", IOobject::NO_READ), np);
236  IOField<scalar> tc(c.fieldIOobject("tc", IOobject::NO_READ), np);
237  IOField<scalar> ms(c.fieldIOobject("ms", IOobject::NO_READ), np);
238  IOField<scalar> injector
239  (
240  c.fieldIOobject("injector", IOobject::NO_READ),
241  np
242  );
243  IOField<scalar> tMom(c.fieldIOobject("tMom", IOobject::NO_READ), np);
244  IOField<scalar> user(c.fieldIOobject("user", IOobject::NO_READ), np);
245 
246  label i = 0;
248  {
249  const SprayParcel<ParcelType>& p = iter();
250  d0[i] = p.d0_;
251  position0[i] = p.position0_;
252  sigma[i] = p.sigma_;
253  mu[i] = p.mu_;
254  liquidCore[i] = p.liquidCore_;
255  KHindex[i] = p.KHindex_;
256  y[i] = p.y_;
257  yDot[i] = p.yDot_;
258  tc[i] = p.tc_;
259  ms[i] = p.ms_;
260  injector[i] = p.injector_;
261  tMom[i] = p.tMom_;
262  user[i] = p.user_;
263  i++;
264  }
265 
266  d0.write();
267  position0.write();
268  sigma.write();
269  mu.write();
270  liquidCore.write();
271  KHindex.write();
272  y.write();
273  yDot.write();
274  tc.write();
275  ms.write();
276  injector.write();
277  tMom.write();
278  user.write();
279 }
280 
281 
282 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
283 
284 template<class ParcelType>
285 Foam::Ostream& Foam::operator<<
286 (
287  Ostream& os,
289 )
290 {
291  if (os.format() == IOstream::ASCII)
292  {
293  os << static_cast<const ParcelType&>(p)
294  << token::SPACE << p.d0()
295  << token::SPACE << p.position0()
296  << token::SPACE << p.sigma()
297  << token::SPACE << p.mu()
298  << token::SPACE << p.liquidCore()
299  << token::SPACE << p.KHindex()
300  << token::SPACE << p.y()
301  << token::SPACE << p.yDot()
302  << token::SPACE << p.tc()
303  << token::SPACE << p.ms()
304  << token::SPACE << p.injector()
305  << token::SPACE << p.tMom()
306  << token::SPACE << p.user();
307  }
308  else
309  {
310  os << static_cast<const ParcelType&>(p);
311  os.write
312  (
313  reinterpret_cast<const char*>(&p.d0_),
315  );
316  }
317 
318  // Check state of Ostream
319  os.check
320  (
321  "Ostream& operator<<(Ostream&, const SprayParcel<ParcelType>&)"
322  );
323 
324  return os;
325 }
326 
327 
328 // ************************************************************************* //
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Foam::SprayParcel::position0_
vector position0_
Injection position.
Definition: SprayParcel.H:139
p
p
Definition: pEqn.H:62
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
Foam::SprayParcel::tc_
scalar tc_
Characteristic time (used in atomization and/or breakup model)
Definition: SprayParcel.H:160
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::SprayParcel::liquidCore_
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet )
Definition: SprayParcel.H:148
Foam::string
A class for handling character strings derived from std::string.
Definition: string.H:74
Foam::SprayParcel::sizeofFields_
static const std::size_t sizeofFields_
Size in bytes of the fields.
Definition: SprayParcel.H:64
Foam::SprayParcel::readFields
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
Definition: SprayParcelIO.C:121
Foam::SprayParcel::d0_
scalar d0_
Initial droplet diameter.
Definition: SprayParcel.H:136
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::SprayParcel::ms_
scalar ms_
Stripped parcel mass due to breakup.
Definition: SprayParcel.H:163
Foam::IOstream::ASCII
@ ASCII
Definition: IOstream.H:88
Foam::SprayParcel::sigma_
scalar sigma_
Liquid surface tension [N/m].
Definition: SprayParcel.H:142
Foam::label
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
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::SprayParcel
Reacing spray parcel, with added functionality for atomization and breakup.
Definition: SprayParcel.H:43
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Foam::SprayParcel::tMom_
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.)
Definition: SprayParcel.H:170
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::SprayParcel::y_
scalar y_
Spherical deviation.
Definition: SprayParcel.H:154
Foam::SprayParcel::KHindex_
scalar KHindex_
Index for KH Breakup.
Definition: SprayParcel.H:151
Foam::SprayParcel::mu_
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:145
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::SprayParcel::user_
scalar user_
Passive scalar (extra variable to be defined by user)
Definition: SprayParcel.H:173
Foam::SprayParcel::writeFields
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
Definition: SprayParcelIO.C:211
Foam::SprayParcel::yDot_
scalar yDot_
Rate of change of spherical deviation.
Definition: SprayParcel.H:157
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
SprayParcel.H
readFields
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
Definition: redistributePar.C:589
Foam::Cloud
Base cloud calls templated on particle type.
Definition: Cloud.H:52
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::SprayParcel::injector_
scalar injector_
Injected from injector (needed e.g. for calculating distance.
Definition: SprayParcel.H:167
Foam::SprayParcel::SprayParcel
SprayParcel(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from owner, position, and cloud owner.
Definition: SprayParcelI.H:108
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::token::SPACE
@ SPACE
Definition: token.H:95