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-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 "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  bool valid = c.size();
90 
92 
93  IOField<vector> f(c.fieldIOobject("f", IOobject::MUST_READ), valid);
94  c.checkFieldIOobject(c, f);
95 
96  IOField<vector> angularMomentum
97  (
98  c.fieldIOobject("angularMomentum", IOobject::MUST_READ),
99  valid
100  );
101  c.checkFieldIOobject(c, angularMomentum);
102 
103  IOField<vector> torque
104  (
105  c.fieldIOobject("torque", IOobject::MUST_READ),
106  valid
107  );
108  c.checkFieldIOobject(c, torque);
109 
110  labelFieldCompactIOField collisionRecordsPairAccessed
111  (
112  c.fieldIOobject("collisionRecordsPairAccessed", IOobject::MUST_READ),
113  valid
114  );
115  c.checkFieldFieldIOobject(c, collisionRecordsPairAccessed);
116 
117  labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
118  (
119  c.fieldIOobject
120  (
121  "collisionRecordsPairOrigProcOfOther",
122  IOobject::MUST_READ
123  ),
124  valid
125  );
126  c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
127 
128  labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
129  (
130  c.fieldIOobject
131  (
132  "collisionRecordsPairOrigIdOfOther",
133  IOobject::MUST_READ
134  ),
135  valid
136  );
137  c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
138 
139  pairDataFieldCompactIOField collisionRecordsPairData
140  (
141  c.fieldIOobject("collisionRecordsPairData", IOobject::MUST_READ),
142  valid
143  );
144  c.checkFieldFieldIOobject(c, collisionRecordsPairData);
145 
146  labelFieldCompactIOField collisionRecordsWallAccessed
147  (
148  c.fieldIOobject("collisionRecordsWallAccessed", IOobject::MUST_READ),
149  valid
150  );
151  c.checkFieldFieldIOobject(c, collisionRecordsWallAccessed);
152 
153  vectorFieldCompactIOField collisionRecordsWallPRel
154  (
155  c.fieldIOobject("collisionRecordsWallPRel", IOobject::MUST_READ),
156  valid
157  );
158  c.checkFieldFieldIOobject(c, collisionRecordsWallPRel);
159 
160  wallDataFieldCompactIOField collisionRecordsWallData
161  (
162  c.fieldIOobject("collisionRecordsWallData", IOobject::MUST_READ),
163  valid
164  );
165  c.checkFieldFieldIOobject(c, collisionRecordsWallData);
166 
167  label i = 0;
168 
169  forAllIter(typename CloudType, c, iter)
170  {
171  CollidingParcel<ParcelType>& p = iter();
172 
173  p.f_ = f[i];
174  p.angularMomentum_ = angularMomentum[i];
175  p.torque_ = torque[i];
176 
178  (
179  collisionRecordsPairAccessed[i],
180  collisionRecordsPairOrigProcOfOther[i],
181  collisionRecordsPairOrigIdOfOther[i],
182  collisionRecordsPairData[i],
183  collisionRecordsWallAccessed[i],
184  collisionRecordsWallPRel[i],
185  collisionRecordsWallData[i]
186  );
187 
188  i++;
189  }
190 }
191 
192 
193 template<class ParcelType>
194 template<class CloudType>
196 {
197  ParcelType::writeFields(c);
198 
199  label np = c.size();
200 
201  IOField<vector> f(c.fieldIOobject("f", IOobject::NO_READ), np);
202  IOField<vector> angularMomentum
203  (
204  c.fieldIOobject("angularMomentum", IOobject::NO_READ),
205  np
206  );
207  IOField<vector> torque(c.fieldIOobject("torque", IOobject::NO_READ), np);
208 
209  labelFieldCompactIOField collisionRecordsPairAccessed
210  (
211  c.fieldIOobject("collisionRecordsPairAccessed", IOobject::NO_READ),
212  np
213  );
214  labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
215  (
216  c.fieldIOobject
217  (
218  "collisionRecordsPairOrigProcOfOther",
219  IOobject::NO_READ
220  ),
221  np
222  );
223  labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
224  (
225  c.fieldIOobject("collisionRecordsPairOrigIdOfOther", IOobject::NO_READ),
226  np
227  );
228  pairDataFieldCompactIOField collisionRecordsPairData
229  (
230  c.fieldIOobject("collisionRecordsPairData", IOobject::NO_READ),
231  np
232  );
233  labelFieldCompactIOField collisionRecordsWallAccessed
234  (
235  c.fieldIOobject("collisionRecordsWallAccessed", IOobject::NO_READ),
236  np
237  );
238  vectorFieldCompactIOField collisionRecordsWallPRel
239  (
240  c.fieldIOobject("collisionRecordsWallPRel", IOobject::NO_READ),
241  np
242  );
243  wallDataFieldCompactIOField collisionRecordsWallData
244  (
245  c.fieldIOobject("collisionRecordsWallData", IOobject::NO_READ),
246  np
247  );
248 
249  label i = 0;
250 
251  forAllConstIter(typename CloudType, c, iter)
252  {
253  const CollidingParcel<ParcelType>& p = iter();
254 
255  f[i] = p.f();
256  angularMomentum[i] = p.angularMomentum();
257  torque[i] = p.torque();
258 
259  collisionRecordsPairAccessed[i] = p.collisionRecords().pairAccessed();
260  collisionRecordsPairOrigProcOfOther[i] =
261  p.collisionRecords().pairOrigProcOfOther();
262  collisionRecordsPairOrigIdOfOther[i] =
263  p.collisionRecords().pairOrigIdOfOther();
264  collisionRecordsPairData[i] = p.collisionRecords().pairData();
265  collisionRecordsWallAccessed[i] = p.collisionRecords().wallAccessed();
266  collisionRecordsWallPRel[i] = p.collisionRecords().wallPRel();
267  collisionRecordsWallData[i] = p.collisionRecords().wallData();
268 
269  i++;
270  }
271 
272  const bool valid = (np > 0);
273 
274  f.write(valid);
275  angularMomentum.write(valid);
276  torque.write(valid);
277 
278  collisionRecordsPairAccessed.write(valid);
279  collisionRecordsPairOrigProcOfOther.write(valid);
280  collisionRecordsPairOrigIdOfOther.write(valid);
281  collisionRecordsPairData.write(valid);
282  collisionRecordsWallAccessed.write(valid);
283  collisionRecordsWallPRel.write(valid);
284  collisionRecordsWallData.write(valid);
285 }
286 
287 
288 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
289 
290 template<class ParcelType>
291 Foam::Ostream& Foam::operator<<
292 (
293  Ostream& os,
295 )
296 {
297  if (os.format() == IOstream::ASCII)
298  {
299  os << static_cast<const ParcelType&>(p)
300  << token::SPACE << p.f_
301  << token::SPACE << p.angularMomentum_
302  << token::SPACE << p.torque_
303  << token::SPACE << p.collisionRecords_;
304  }
305  else
306  {
307  os << static_cast<const ParcelType&>(p);
308  os.write
309  (
310  reinterpret_cast<const char*>(&p.f_),
312  );
313  os << p.collisionRecords();
314  }
315 
316  // Check state of Ostream
317  os.check
318  (
319  "Ostream& operator<<(Ostream&, const CollidingParcel<ParcelType>&)"
320  );
321 
322  return os;
323 }
324 
325 
326 // ************************************************************************* //
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
CollisionRecordList< vector, vector > collisionRecordList
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
Wrapper around kinematic parcel types to add collision modelling.
#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
void checkFieldFieldIOobject(const Cloud< ParticleType > &c, const CompactIOField< Field< DataType >, DataType > &data) const
Check lagrangian data fieldfield.
Definition: CloudIO.C:205
const vector & f() const
Return const access to force.
static void writeFields(const CloudType &c)
Write.
label size() const
Definition: Cloud.H:162
collisionRecordList collisionRecords_
Particle collision records.
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
virtual Istream & read(token &)=0
Return next token from stream.
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
const vector & torque() const
Return const access to torque.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const vector & angularMomentum() const
Return const access to angular momentum.
labelList f(nPoints)
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
vector torque_
Torque on particle due to collisions in global.
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.
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
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
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
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
static void readFields(CloudType &c)
Read.