CollidingParcelIO.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 "CollidingParcel.H"
27 #include "IOstreams.H"
28 #include "IOField.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 template<class ParcelType>
33 Foam::string Foam::CollidingParcel<ParcelType>::propertyList_ =
35 
36 template<class ParcelType>
38 (
39  offsetof(CollidingParcel<ParcelType>, collisionRecords_)
40  - offsetof(CollidingParcel<ParcelType>, f_)
41 );
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class ParcelType>
48 (
49  const polyMesh& mesh,
50  Istream& is,
51  bool readFields
52 )
53 :
54  ParcelType(mesh, is, readFields),
55  f_(Zero),
56  angularMomentum_(Zero),
57  torque_(Zero),
58  collisionRecords_()
59 {
60  if (readFields)
61  {
62  if (is.format() == IOstream::ASCII)
63  {
64  is >> f_;
65  is >> angularMomentum_;
66  is >> torque_;
67  }
68  else
69  {
70  is.read(reinterpret_cast<char*>(&f_), sizeofFields_);
71  }
72 
73  is >> collisionRecords_;
74  }
75 
76  // Check state of Istream
77  is.check
78  (
79  "CollidingParcel<ParcelType>::Collisions"
80  "(const polyMesh&, Istream&, bool)"
81  );
82 }
83 
84 
85 template<class ParcelType>
86 template<class CloudType>
88 {
89  if (!c.size())
90  {
91  return;
92  }
93 
95 
96  IOField<vector> f(c.fieldIOobject("f", IOobject::MUST_READ));
97  c.checkFieldIOobject(c, f);
98 
99  IOField<vector> angularMomentum
100  (
101  c.fieldIOobject("angularMomentum", IOobject::MUST_READ)
102  );
103  c.checkFieldIOobject(c, angularMomentum);
104 
105  IOField<vector> torque(c.fieldIOobject("torque", IOobject::MUST_READ));
106  c.checkFieldIOobject(c, torque);
107 
108  labelFieldCompactIOField collisionRecordsPairAccessed
109  (
110  c.fieldIOobject("collisionRecordsPairAccessed", IOobject::MUST_READ)
111  );
112  c.checkFieldFieldIOobject(c, collisionRecordsPairAccessed);
113 
114  labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
115  (
116  c.fieldIOobject
117  (
118  "collisionRecordsPairOrigProcOfOther",
119  IOobject::MUST_READ
120  )
121  );
122  c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
123 
124  labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
125  (
126  c.fieldIOobject
127  (
128  "collisionRecordsPairOrigIdOfOther",
129  IOobject::MUST_READ
130  )
131  );
132  c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
133 
134  pairDataFieldCompactIOField collisionRecordsPairData
135  (
136  c.fieldIOobject("collisionRecordsPairData", IOobject::MUST_READ)
137  );
138  c.checkFieldFieldIOobject(c, collisionRecordsPairData);
139 
140  labelFieldCompactIOField collisionRecordsWallAccessed
141  (
142  c.fieldIOobject("collisionRecordsWallAccessed", IOobject::MUST_READ)
143  );
144  c.checkFieldFieldIOobject(c, collisionRecordsWallAccessed);
145 
146  vectorFieldCompactIOField collisionRecordsWallPRel
147  (
148  c.fieldIOobject("collisionRecordsWallPRel", IOobject::MUST_READ)
149  );
150  c.checkFieldFieldIOobject(c, collisionRecordsWallPRel);
151 
152  wallDataFieldCompactIOField collisionRecordsWallData
153  (
154  c.fieldIOobject("collisionRecordsWallData", IOobject::MUST_READ)
155  );
156  c.checkFieldFieldIOobject(c, collisionRecordsWallData);
157 
158  label i = 0;
159 
160  forAllIter(typename CloudType, c, iter)
161  {
162  CollidingParcel<ParcelType>& p = iter();
163 
164  p.f_ = f[i];
165  p.angularMomentum_ = angularMomentum[i];
166  p.torque_ = torque[i];
167 
169  (
170  collisionRecordsPairAccessed[i],
171  collisionRecordsPairOrigProcOfOther[i],
172  collisionRecordsPairOrigIdOfOther[i],
173  collisionRecordsPairData[i],
174  collisionRecordsWallAccessed[i],
175  collisionRecordsWallPRel[i],
176  collisionRecordsWallData[i]
177  );
178 
179  i++;
180  }
181 }
182 
183 
184 template<class ParcelType>
185 template<class CloudType>
187 {
188  ParcelType::writeFields(c);
189 
190  label np = c.size();
191 
192  IOField<vector> f(c.fieldIOobject("f", IOobject::NO_READ), np);
193  IOField<vector> angularMomentum
194  (
195  c.fieldIOobject("angularMomentum", IOobject::NO_READ),
196  np
197  );
198  IOField<vector> torque(c.fieldIOobject("torque", IOobject::NO_READ), np);
199 
200  labelFieldCompactIOField collisionRecordsPairAccessed
201  (
202  c.fieldIOobject("collisionRecordsPairAccessed", IOobject::NO_READ),
203  np
204  );
205  labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
206  (
207  c.fieldIOobject
208  (
209  "collisionRecordsPairOrigProcOfOther",
210  IOobject::NO_READ
211  ),
212  np
213  );
214  labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
215  (
216  c.fieldIOobject("collisionRecordsPairOrigIdOfOther", IOobject::NO_READ),
217  np
218  );
219  pairDataFieldCompactIOField collisionRecordsPairData
220  (
221  c.fieldIOobject("collisionRecordsPairData", IOobject::NO_READ),
222  np
223  );
224  labelFieldCompactIOField collisionRecordsWallAccessed
225  (
226  c.fieldIOobject("collisionRecordsWallAccessed", IOobject::NO_READ),
227  np
228  );
229  vectorFieldCompactIOField collisionRecordsWallPRel
230  (
231  c.fieldIOobject("collisionRecordsWallPRel", IOobject::NO_READ),
232  np
233  );
234  wallDataFieldCompactIOField collisionRecordsWallData
235  (
236  c.fieldIOobject("collisionRecordsWallData", IOobject::NO_READ),
237  np
238  );
239 
240  label i = 0;
241 
242  forAllConstIter(typename CloudType, c, iter)
243  {
244  const CollidingParcel<ParcelType>& p = iter();
245 
246  f[i] = p.f();
247  angularMomentum[i] = p.angularMomentum();
248  torque[i] = p.torque();
249 
250  collisionRecordsPairAccessed[i] = p.collisionRecords().pairAccessed();
251  collisionRecordsPairOrigProcOfOther[i] =
252  p.collisionRecords().pairOrigProcOfOther();
253  collisionRecordsPairOrigIdOfOther[i] =
254  p.collisionRecords().pairOrigIdOfOther();
255  collisionRecordsPairData[i] = p.collisionRecords().pairData();
256  collisionRecordsWallAccessed[i] = p.collisionRecords().wallAccessed();
257  collisionRecordsWallPRel[i] = p.collisionRecords().wallPRel();
258  collisionRecordsWallData[i] = p.collisionRecords().wallData();
259 
260  i++;
261  }
262 
263  f.write();
264  angularMomentum.write();
265  torque.write();
266 
267  collisionRecordsPairAccessed.write();
268  collisionRecordsPairOrigProcOfOther.write();
269  collisionRecordsPairOrigIdOfOther.write();
270  collisionRecordsPairData.write();
271  collisionRecordsWallAccessed.write();
272  collisionRecordsWallPRel.write();
273  collisionRecordsWallData.write();
274 }
275 
276 
277 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
278 
279 template<class ParcelType>
280 Foam::Ostream& Foam::operator<<
281 (
282  Ostream& os,
284 )
285 {
286  if (os.format() == IOstream::ASCII)
287  {
288  os << static_cast<const ParcelType&>(p)
289  << token::SPACE << p.f_
290  << token::SPACE << p.angularMomentum_
291  << token::SPACE << p.torque_
292  << token::SPACE << p.collisionRecords_;
293  }
294  else
295  {
296  os << static_cast<const ParcelType&>(p);
297  os.write
298  (
299  reinterpret_cast<const char*>(&p.f_),
301  );
302  os << p.collisionRecords();
303  }
304 
305  // Check state of Ostream
306  os.check
307  (
308  "Ostream& operator<<(Ostream&, const CollidingParcel<ParcelType>&)"
309  );
310 
311  return os;
312 }
313 
314 
315 // ************************************************************************* //
const vector & torque() const
Return const access to torque.
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
vector f_
Force on particle due to collisions [N].
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 checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:215
CollisionRecordList< vector, vector > collisionRecordList
Wrapper around kinematic parcel types to add collision modelling.
void checkFieldFieldIOobject(const Cloud< ParticleType > &c, const CompactIOField< Field< DataType >, DataType > &data) const
Check lagrangian data fieldfield.
Definition: CloudIO.C:234
#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
static void writeFields(const CloudType &c)
Write.
const vector & f() const
Return const access to force.
collisionRecordList collisionRecords_
Particle collision records.
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...
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
virtual Istream & read(token &)=0
Return next token from stream.
static const zero Zero
Definition: zero.H:91
CollidingParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
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:53
labelList f(nPoints)
label size() const
Definition: Cloud.H:175
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
vector torque_
Torque on particle due to collisions in global.
const vector & angularMomentum() const
Return const access to angular momentum.
A Field of objects of type <T> with automated input and output using a compact storage. Behaves like IOField except when binary output in case it writes a CompactListList.
const dimensionedScalar c
Speed of light in a vacuum.
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
static void readFields(CloudType &c)
Read.