KinematicParcelIO.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-2017 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 "KinematicParcel.H"
27 #include "IOstreams.H"
28 #include "IOField.H"
29 #include "Cloud.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 template<class ParcelType>
34 Foam::string Foam::KinematicParcel<ParcelType>::propertyList_ =
36 
37 template<class ParcelType>
39 (
40  offsetof(KinematicParcel<ParcelType>, rhoc_)
41  - offsetof(KinematicParcel<ParcelType>, active_)
42 );
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 template<class ParcelType>
49 (
50  const polyMesh& mesh,
51  Istream& is,
52  bool readFields
53 )
54 :
55  ParcelType(mesh, is, readFields),
56  active_(false),
57  typeId_(0),
58  nParticle_(0.0),
59  d_(0.0),
60  dTarget_(0.0),
61  U_(Zero),
62  rho_(0.0),
63  age_(0.0),
64  tTurb_(0.0),
65  UTurb_(Zero),
66  rhoc_(0.0),
67  Uc_(Zero),
68  muc_(0.0)
69 {
70  if (readFields)
71  {
72  if (is.format() == IOstream::ASCII)
73  {
74  active_ = readBool(is);
75  typeId_ = readLabel(is);
76  nParticle_ = readScalar(is);
77  d_ = readScalar(is);
78  dTarget_ = readScalar(is);
79  is >> U_;
80  rho_ = readScalar(is);
81  age_ = readScalar(is);
82  tTurb_ = readScalar(is);
83  is >> UTurb_;
84  }
85  else
86  {
87  is.read(reinterpret_cast<char*>(&active_), sizeofFields_);
88  }
89  }
90 
91  // Check state of Istream
92  is.check
93  (
94  "KinematicParcel<ParcelType>::KinematicParcel"
95  "(const polyMesh&, Istream&, bool)"
96  );
97 }
98 
99 
100 template<class ParcelType>
101 template<class CloudType>
103 {
104  bool valid = c.size();
105 
107 
108  IOField<label> active
109  (
110  c.fieldIOobject("active", IOobject::MUST_READ),
111  valid
112  );
113  c.checkFieldIOobject(c, active);
114 
115  IOField<label> typeId
116  (
117  c.fieldIOobject("typeId", IOobject::MUST_READ),
118  valid
119  );
120  c.checkFieldIOobject(c, typeId);
121 
122  IOField<scalar> nParticle
123  (
124  c.fieldIOobject("nParticle", IOobject::MUST_READ),
125  valid
126  );
127  c.checkFieldIOobject(c, nParticle);
128 
130  (
131  c.fieldIOobject("d", IOobject::MUST_READ),
132  valid
133  );
134  c.checkFieldIOobject(c, d);
135 
136  IOField<scalar> dTarget
137  (
138  c.fieldIOobject("dTarget", IOobject::MUST_READ),
139  valid
140  );
141  c.checkFieldIOobject(c, dTarget);
142 
144  (
145  c.fieldIOobject("U", IOobject::MUST_READ),
146  valid
147  );
148  c.checkFieldIOobject(c, U);
149 
151  (
152  c.fieldIOobject("rho", IOobject::MUST_READ),
153  valid
154  );
155  c.checkFieldIOobject(c, rho);
156 
157  IOField<scalar> age
158  (
159  c.fieldIOobject("age", IOobject::MUST_READ),
160  valid
161  );
162  c.checkFieldIOobject(c, age);
163 
164  IOField<scalar> tTurb
165  (
166  c.fieldIOobject("tTurb", IOobject::MUST_READ),
167  valid
168  );
169  c.checkFieldIOobject(c, tTurb);
170 
171  IOField<vector> UTurb
172  (
173  c.fieldIOobject("UTurb", IOobject::MUST_READ),
174  valid
175  );
176  c.checkFieldIOobject(c, UTurb);
177 
178  label i = 0;
179 
180  forAllIter(typename CloudType, c, iter)
181  {
182  KinematicParcel<ParcelType>& p = iter();
183 
184  p.active_ = active[i];
185  p.typeId_ = typeId[i];
186  p.nParticle_ = nParticle[i];
187  p.d_ = d[i];
188  p.dTarget_ = dTarget[i];
189  p.U_ = U[i];
190  p.rho_ = rho[i];
191  p.age_ = age[i];
192  p.tTurb_ = tTurb[i];
193  p.UTurb_ = UTurb[i];
194 
195  i++;
196  }
197 }
198 
199 
200 template<class ParcelType>
201 template<class CloudType>
203 {
204  ParcelType::writeFields(c);
205 
206  label np = c.size();
207 
208  IOField<label> active(c.fieldIOobject("active", IOobject::NO_READ), np);
209  IOField<label> typeId(c.fieldIOobject("typeId", IOobject::NO_READ), np);
210  IOField<scalar> nParticle
211  (
212  c.fieldIOobject("nParticle", IOobject::NO_READ),
213  np
214  );
215  IOField<scalar> d(c.fieldIOobject("d", IOobject::NO_READ), np);
216  IOField<scalar> dTarget(c.fieldIOobject("dTarget", IOobject::NO_READ), np);
217  IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
218  IOField<scalar> rho(c.fieldIOobject("rho", IOobject::NO_READ), np);
219  IOField<scalar> age(c.fieldIOobject("age", IOobject::NO_READ), np);
220  IOField<scalar> tTurb(c.fieldIOobject("tTurb", IOobject::NO_READ), np);
221  IOField<vector> UTurb(c.fieldIOobject("UTurb", IOobject::NO_READ), np);
222 
223  label i = 0;
224 
225  forAllConstIter(typename CloudType, c, iter)
226  {
227  const KinematicParcel<ParcelType>& p = iter();
228 
229  active[i] = p.active();
230  typeId[i] = p.typeId();
231  nParticle[i] = p.nParticle();
232  d[i] = p.d();
233  dTarget[i] = p.dTarget();
234  U[i] = p.U();
235  rho[i] = p.rho();
236  age[i] = p.age();
237  tTurb[i] = p.tTurb();
238  UTurb[i] = p.UTurb();
239 
240  i++;
241  }
242 
243  const bool valid = np > 0;
244 
245  active.write(valid);
246  typeId.write(valid);
247  nParticle.write(valid);
248  d.write(valid);
249  dTarget.write(valid);
250  U.write(valid);
251  rho.write(valid);
252  age.write(valid);
253  tTurb.write(valid);
254  UTurb.write(valid);
255 }
256 
257 
258 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
259 
260 template<class ParcelType>
261 Foam::Ostream& Foam::operator<<
262 (
263  Ostream& os,
265 )
266 {
267  if (os.format() == IOstream::ASCII)
268  {
269  os << static_cast<const ParcelType&>(p)
270  << token::SPACE << p.active()
271  << token::SPACE << p.typeId()
272  << token::SPACE << p.nParticle()
273  << token::SPACE << p.d()
274  << token::SPACE << p.dTarget()
275  << token::SPACE << p.U()
276  << token::SPACE << p.rho()
277  << token::SPACE << p.age()
278  << token::SPACE << p.tTurb()
279  << token::SPACE << p.UTurb();
280  }
281  else
282  {
283  os << static_cast<const ParcelType&>(p);
284  os.write
285  (
286  reinterpret_cast<const char*>(&p.active_),
288  );
289  }
290 
291  // Check state of Ostream
292  os.check
293  (
294  "Ostream& operator<<(Ostream&, const KinematicParcel<ParcelType>&)"
295  );
296 
297  return os;
298 }
299 
300 
301 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
const vector & U() const
Return const access to velocity.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
U
Definition: pEqn.H:83
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
static void readFields(CloudType &c)
Read.
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
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
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
scalar nParticle() const
Return const access to number of particles.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar d() const
Return const access to diameter.
scalar tTurb_
Time spent in turbulent eddy [s].
scalar age_
Age [s].
scalar rho() const
Return const access to density.
label size() const
Definition: Cloud.H:162
bool readBool(Istream &)
Definition: boolIO.C:60
KinematicParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar dTarget() const
Return const access to target diameter.
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
scalar dTarget_
Target diameter [m].
scalar age() const
Return const access to the age.
virtual Istream & read(token &)=0
Return next token from stream.
scalar rho_
Density [kg/m3].
bool active_
Active flag - tracking inactive when active = false.
label typeId_
Parcel type id.
vector UTurb_
Turbulent velocity fluctuation [m/s].
static const zero Zero
Definition: zero.H:91
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
bool active() const
Return const access to active flag.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
scalar nParticle_
Number of particles in Parcel.
vector U_
Velocity of Parcel [m/s].
const dimensionedScalar c
Speed of light in a vacuum.
static void writeFields(const CloudType &c)
Write.
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:186
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
virtual bool write(const bool valid=true) const
Write using setting from DB.
label typeId() const
Return const access to type id.
A class for handling character strings derived from std::string.
Definition: string.H:74
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:166
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
scalar d_
Diameter [m].