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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "CollidingParcel.H"
30#include "IOstreams.H"
31#include "IOField.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35template<class ParcelType>
38
39template<class ParcelType>
41(
42 offsetof(CollidingParcel<ParcelType>, collisionRecords_)
43 - offsetof(CollidingParcel<ParcelType>, f_)
44);
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49template<class ParcelType>
51(
52 const polyMesh& mesh,
53 Istream& is,
54 bool readFields,
55 bool newFormat
56)
57:
58 ParcelType(mesh, is, readFields, newFormat),
59 f_(Zero),
60 angularMomentum_(Zero),
61 torque_(Zero),
62 collisionRecords_()
63{
64 if (readFields)
65 {
66 if (is.format() == IOstream::ASCII)
67 {
68 is >> f_;
69 is >> angularMomentum_;
70 is >> torque_;
71 }
72 else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
73 {
74 // Non-native label or scalar size
75
76 is.beginRawRead();
77
78 readRawScalar(is, f_.data(), vector::nComponents);
79 readRawScalar(is, angularMomentum_.data(), vector::nComponents);
80 readRawScalar(is, torque_.data(), vector::nComponents);
81
82 is.endRawRead();
83 }
84 else
85 {
86 is.read(reinterpret_cast<char*>(&f_), sizeofFields);
87 }
88
90 }
91
93}
94
95
96template<class ParcelType>
97template<class CloudType>
99{
100 const bool valid = c.size();
101
103
104 IOField<vector> f(c.fieldIOobject("f", IOobject::MUST_READ), valid);
105 c.checkFieldIOobject(c, f);
106
107 IOField<vector> angularMomentum
108 (
109 c.fieldIOobject("angularMomentum", IOobject::MUST_READ),
110 valid
111 );
112 c.checkFieldIOobject(c, angularMomentum);
113
114 IOField<vector> torque
115 (
116 c.fieldIOobject("torque", IOobject::MUST_READ),
117 valid
118 );
119 c.checkFieldIOobject(c, torque);
120
121 labelFieldCompactIOField collisionRecordsPairAccessed
122 (
123 c.fieldIOobject("collisionRecordsPairAccessed", IOobject::MUST_READ),
124 valid
125 );
126 c.checkFieldFieldIOobject(c, collisionRecordsPairAccessed);
127
128 labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
129 (
130 c.fieldIOobject
131 (
132 "collisionRecordsPairOrigProcOfOther",
134 ),
135 valid
136 );
137 c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
138
139 labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
140 (
141 c.fieldIOobject
142 (
143 "collisionRecordsPairOrigIdOfOther",
145 ),
146 valid
147 );
148 c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
149
150 pairDataFieldCompactIOField collisionRecordsPairData
151 (
152 c.fieldIOobject("collisionRecordsPairData", IOobject::MUST_READ),
153 valid
154 );
155 c.checkFieldFieldIOobject(c, collisionRecordsPairData);
156
157 labelFieldCompactIOField collisionRecordsWallAccessed
158 (
159 c.fieldIOobject("collisionRecordsWallAccessed", IOobject::MUST_READ),
160 valid
161 );
162 c.checkFieldFieldIOobject(c, collisionRecordsWallAccessed);
163
164 vectorFieldCompactIOField collisionRecordsWallPRel
165 (
166 c.fieldIOobject("collisionRecordsWallPRel", IOobject::MUST_READ),
167 valid
168 );
169 c.checkFieldFieldIOobject(c, collisionRecordsWallPRel);
170
171 wallDataFieldCompactIOField collisionRecordsWallData
172 (
173 c.fieldIOobject("collisionRecordsWallData", IOobject::MUST_READ),
174 valid
175 );
176 c.checkFieldFieldIOobject(c, collisionRecordsWallData);
177
178 label i = 0;
179
181 {
182 p.f_ = f[i];
183 p.angularMomentum_ = angularMomentum[i];
184 p.torque_ = torque[i];
185
186 p.collisionRecords_ = collisionRecordList
187 (
188 collisionRecordsPairAccessed[i],
189 collisionRecordsPairOrigProcOfOther[i],
190 collisionRecordsPairOrigIdOfOther[i],
191 collisionRecordsPairData[i],
192 collisionRecordsWallAccessed[i],
193 collisionRecordsWallPRel[i],
194 collisionRecordsWallData[i]
195 );
196
197 ++i;
198 }
199}
200
201
202template<class ParcelType>
203template<class CloudType>
205{
207
208 const label np = c.size();
209 const bool valid = np;
210
211
212 IOField<vector> f(c.fieldIOobject("f", IOobject::NO_READ), np);
213 IOField<vector> angMom
214 (
215 c.fieldIOobject("angularMomentum", IOobject::NO_READ),
216 np
217 );
218 IOField<vector> torque(c.fieldIOobject("torque", IOobject::NO_READ), np);
219
220 labelFieldCompactIOField collisionRecordsPairAccessed
221 (
222 c.fieldIOobject("collisionRecordsPairAccessed", IOobject::NO_READ),
223 np
224 );
225 labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
226 (
227 c.fieldIOobject
228 (
229 "collisionRecordsPairOrigProcOfOther",
231 ),
232 np
233 );
234 labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
235 (
236 c.fieldIOobject("collisionRecordsPairOrigIdOfOther", IOobject::NO_READ),
237 np
238 );
239 pairDataFieldCompactIOField collisionRecordsPairData
240 (
241 c.fieldIOobject("collisionRecordsPairData", IOobject::NO_READ),
242 np
243 );
244 labelFieldCompactIOField collisionRecordsWallAccessed
245 (
246 c.fieldIOobject("collisionRecordsWallAccessed", IOobject::NO_READ),
247 np
248 );
249 vectorFieldCompactIOField collisionRecordsWallPRel
250 (
251 c.fieldIOobject("collisionRecordsWallPRel", IOobject::NO_READ),
252 np
253 );
254 wallDataFieldCompactIOField collisionRecordsWallData
255 (
256 c.fieldIOobject("collisionRecordsWallData", IOobject::NO_READ),
257 np
258 );
259
260 label i = 0;
261 for (const CollidingParcel<ParcelType>& p : c)
262 {
263 f[i] = p.f();
264 angMom[i] = p.angularMomentum();
265 torque[i] = p.torque();
266
267 collisionRecordsPairAccessed[i] = p.collisionRecords().pairAccessed();
268 collisionRecordsPairOrigProcOfOther[i] =
269 p.collisionRecords().pairOrigProcOfOther();
270 collisionRecordsPairOrigIdOfOther[i] =
271 p.collisionRecords().pairOrigIdOfOther();
272 collisionRecordsPairData[i] = p.collisionRecords().pairData();
273 collisionRecordsWallAccessed[i] = p.collisionRecords().wallAccessed();
274 collisionRecordsWallPRel[i] = p.collisionRecords().wallPRel();
275 collisionRecordsWallData[i] = p.collisionRecords().wallData();
276
277 ++i;
278 }
279
280 f.write(valid);
281 angMom.write(valid);
282 torque.write(valid);
283
284 collisionRecordsPairAccessed.write(valid);
285 collisionRecordsPairOrigProcOfOther.write(valid);
286 collisionRecordsPairOrigIdOfOther.write(valid);
287 collisionRecordsPairData.write(valid);
288 collisionRecordsWallAccessed.write(valid);
289 collisionRecordsWallPRel.write(valid);
290 collisionRecordsWallData.write(valid);
291}
292
293
294template<class ParcelType>
296(
297 Ostream& os,
298 const wordRes& filters,
299 const word& delim,
300 const bool namesOnly
301) const
302{
303 ParcelType::writeProperties(os, filters, delim, namesOnly);
304
305 #undef writeProp
306 #define writeProp(Name, Value) \
307 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
308
309 writeProp("f", f_);
310 writeProp("angularMomentum", angularMomentum_);
311 writeProp("torque", torque_);
312 //writeProp("collisionRecords", collisionRecords_);
313
314 #undef writeProp
315}
316
317
318template<class ParcelType>
319template<class CloudType>
321(
322 CloudType& c,
323 const objectRegistry& obr
324)
325{
327
328 if (!c.size()) return;
329
330 const auto& f = cloud::lookupIOField<vector>("f", obr);
331 const auto& angMom = cloud::lookupIOField<vector>("angularMomentum", obr);
332 const auto& torque = cloud::lookupIOField<vector>("torque", obr);
333
334 label i = 0;
336 {
337 p.f_ = f[i];
338 p.angularMomentum_ = angMom[i];
339 p.torque_ = torque[i];
340
341 ++i;
342 }
343}
344
345
346template<class ParcelType>
347template<class CloudType>
349(
350 const CloudType& c,
351 objectRegistry& obr
352)
353{
355
356 const label np = c.size();
357
358 auto& f = cloud::createIOField<vector>("f", np, obr);
359 auto& angMom = cloud::createIOField<vector>("angularMomentum", np, obr);
360 auto& torque = cloud::createIOField<vector>("torque", np, obr);
361
362 label i = 0;
363 for (const CollidingParcel<ParcelType>& p : c)
364 {
365 f[i] = p.f();
366 angMom[i] = p.angularMomentum();
367 torque[i] = p.torque();
368
369 ++i;
370 }
371}
372
373
374// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
375
376template<class ParcelType>
377Foam::Ostream& Foam::operator<<
378(
379 Ostream& os,
381)
382{
383 if (os.format() == IOstream::ASCII)
384 {
385 os << static_cast<const ParcelType&>(p)
386 << token::SPACE << p.f_
387 << token::SPACE << p.angularMomentum_
388 << token::SPACE << p.torque_
389 << token::SPACE << p.collisionRecords_;
390 }
391 else
392 {
393 os << static_cast<const ParcelType&>(p);
394 os.write
395 (
396 reinterpret_cast<const char*>(&p.f_),
398 );
399 os << p.collisionRecords();
400 }
401
403 return os;
404}
405
406
407// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Wrapper around kinematic parcel types to add collision modelling.
static const std::size_t sizeofFields
Size in bytes of the fields.
vector f_
Force on particle due to collisions [N].
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
collisionRecordList collisionRecords_
Particle collision records.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
static void readFields(CloudType &c)
Read.
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....
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
streamFormat format() const noexcept
Get the current stream format.
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const noexcept
Definition: IOstream.H:300
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const noexcept
Definition: IOstream.H:309
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
virtual bool endRawRead()=0
End of low-level raw binary read.
virtual bool beginRawRead()=0
Start of low-level raw binary read.
virtual Istream & read(token &)=0
Return next token from stream.
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
Allows specification of different writing frequency of objects registered to the database.
Definition: writeObjects.H:142
static void readObjects(Cloud< injectedParticle > &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
Registry of regIOobjects.
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
static string propertyList()
Definition: particle.H:372
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for handling character strings derived from std::string.
Definition: string.H:79
@ SPACE
Space [isspace].
Definition: token.H:125
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define FUNCTION_NAME
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
CollisionRecordList< vector, vector > collisionRecordList
labelList f(nPoints)