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-2016 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  if (!c.size())
105  {
106  return;
107  }
108 
110 
111  IOField<label> active(c.fieldIOobject("active", IOobject::MUST_READ));
112  c.checkFieldIOobject(c, active);
113 
114  IOField<label> typeId(c.fieldIOobject("typeId", IOobject::MUST_READ));
115  c.checkFieldIOobject(c, typeId);
116 
118  nParticle(c.fieldIOobject("nParticle", IOobject::MUST_READ));
119  c.checkFieldIOobject(c, nParticle);
120 
121  IOField<scalar> d(c.fieldIOobject("d", IOobject::MUST_READ));
122  c.checkFieldIOobject(c, d);
123 
124  IOField<scalar> dTarget(c.fieldIOobject("dTarget", IOobject::MUST_READ));
125  c.checkFieldIOobject(c, dTarget);
126 
127  IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ));
128  c.checkFieldIOobject(c, U);
129 
130  IOField<scalar> rho(c.fieldIOobject("rho", IOobject::MUST_READ));
131  c.checkFieldIOobject(c, rho);
132 
133  IOField<scalar> age(c.fieldIOobject("age", IOobject::MUST_READ));
134  c.checkFieldIOobject(c, age);
135 
136  IOField<scalar> tTurb(c.fieldIOobject("tTurb", IOobject::MUST_READ));
137  c.checkFieldIOobject(c, tTurb);
138 
139  IOField<vector> UTurb(c.fieldIOobject("UTurb", IOobject::MUST_READ));
140  c.checkFieldIOobject(c, UTurb);
141 
142  label i = 0;
143 
144  forAllIter(typename CloudType, c, iter)
145  {
146  KinematicParcel<ParcelType>& p = iter();
147 
148  p.active_ = active[i];
149  p.typeId_ = typeId[i];
150  p.nParticle_ = nParticle[i];
151  p.d_ = d[i];
152  p.dTarget_ = dTarget[i];
153  p.U_ = U[i];
154  p.rho_ = rho[i];
155  p.age_ = age[i];
156  p.tTurb_ = tTurb[i];
157  p.UTurb_ = UTurb[i];
158 
159  i++;
160  }
161 }
162 
163 
164 template<class ParcelType>
165 template<class CloudType>
167 {
168  ParcelType::writeFields(c);
169 
170  label np = c.size();
171 
172  IOField<label> active(c.fieldIOobject("active", IOobject::NO_READ), np);
173  IOField<label> typeId(c.fieldIOobject("typeId", IOobject::NO_READ), np);
174  IOField<scalar> nParticle
175  (
176  c.fieldIOobject("nParticle", IOobject::NO_READ),
177  np
178  );
179  IOField<scalar> d(c.fieldIOobject("d", IOobject::NO_READ), np);
180  IOField<scalar> dTarget(c.fieldIOobject("dTarget", IOobject::NO_READ), np);
181  IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
182  IOField<scalar> rho(c.fieldIOobject("rho", IOobject::NO_READ), np);
183  IOField<scalar> age(c.fieldIOobject("age", IOobject::NO_READ), np);
184  IOField<scalar> tTurb(c.fieldIOobject("tTurb", IOobject::NO_READ), np);
185  IOField<vector> UTurb(c.fieldIOobject("UTurb", IOobject::NO_READ), np);
186 
187  label i = 0;
188 
189  forAllConstIter(typename CloudType, c, iter)
190  {
191  const KinematicParcel<ParcelType>& p = iter();
192 
193  active[i] = p.active();
194  typeId[i] = p.typeId();
195  nParticle[i] = p.nParticle();
196  d[i] = p.d();
197  dTarget[i] = p.dTarget();
198  U[i] = p.U();
199  rho[i] = p.rho();
200  age[i] = p.age();
201  tTurb[i] = p.tTurb();
202  UTurb[i] = p.UTurb();
203 
204  i++;
205  }
206 
207  active.write();
208  typeId.write();
209  nParticle.write();
210  d.write();
211  dTarget.write();
212  U.write();
213  rho.write();
214  age.write();
215  tTurb.write();
216  UTurb.write();
217 }
218 
219 
220 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
221 
222 template<class ParcelType>
223 Foam::Ostream& Foam::operator<<
224 (
225  Ostream& os,
227 )
228 {
229  if (os.format() == IOstream::ASCII)
230  {
231  os << static_cast<const ParcelType&>(p)
232  << token::SPACE << p.active()
233  << token::SPACE << p.typeId()
234  << token::SPACE << p.nParticle()
235  << token::SPACE << p.d()
236  << token::SPACE << p.dTarget()
237  << token::SPACE << p.U()
238  << token::SPACE << p.rho()
239  << token::SPACE << p.age()
240  << token::SPACE << p.tTurb()
241  << token::SPACE << p.UTurb();
242  }
243  else
244  {
245  os << static_cast<const ParcelType&>(p);
246  os.write
247  (
248  reinterpret_cast<const char*>(&p.active_),
250  );
251  }
252 
253  // Check state of Ostream
254  os.check
255  (
256  "Ostream& operator<<(Ostream&, const KinematicParcel<ParcelType>&)"
257  );
258 
259  return os;
260 }
261 
262 
263 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
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.
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:215
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
#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
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
scalar rho() const
Return const access to density.
scalar tTurb_
Time spent in turbulent eddy [s].
scalar age_
Age [s].
KinematicParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
bool readBool(Istream &)
Definition: boolIO.C:60
scalar dTarget() const
Return const access to target diameter.
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
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
scalar dTarget_
Target diameter [m].
scalar tTurb() const
Return const access to time spent in turbulent eddy.
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
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
scalar age() const
Return const access to the age.
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.
label size() const
Definition: Cloud.H:175
vector U_
Velocity of Parcel [m/s].
bool active() const
Return const access to active flag.
const dimensionedScalar c
Speed of light in a vacuum.
virtual bool write() const
Write using setting from DB.
static void writeFields(const CloudType &c)
Write.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:195
volScalarField & p
A class for handling character strings derived from std::string.
Definition: string.H:74
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
scalar d() const
Return const access to diameter.
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
scalar d_
Diameter [m].
scalar nParticle() const
Return const access to number of particles.
const vector & U() const
Return const access to velocity.
label typeId() const
Return const access to type id.