Loading...
Searching...
No Matches
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-2022 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_)
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),
63{
64 if (readFields)
65 {
66 if (is.format() == IOstreamOption::ASCII)
67 {
68 is >> f_;
69 is >> angularMomentum_;
70 is >> torque_;
71 }
72 else
73 {
74 // No non-native streaming
76
77 is.read(reinterpret_cast<char*>(&f_), sizeofFields);
78 }
79
81 }
84}
85
86
87template<class ParcelType>
88template<class CloudType>
90{
91 const bool readOnProc = c.size();
92
93 ParcelType::readFields(c);
94
95 IOField<vector> f(c.newIOobject("f", IOobject::MUST_READ), readOnProc);
96 c.checkFieldIOobject(c, f);
97
98 IOField<vector> angularMomentum
99 (
100 c.newIOobject("angularMomentum", IOobject::MUST_READ),
101 readOnProc
102 );
103 c.checkFieldIOobject(c, angularMomentum);
104
105 IOField<vector> torque
106 (
107 c.newIOobject("torque", IOobject::MUST_READ),
108 readOnProc
109 );
110 c.checkFieldIOobject(c, torque);
111
112 labelFieldCompactIOField collisionRecordsPairAccessed
113 (
114 c.newIOobject("collisionRecordsPairAccessed", IOobject::MUST_READ),
115 readOnProc
116 );
117 c.checkFieldFieldIOobject(c, collisionRecordsPairAccessed);
118
119 labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
120 (
121 c.newIOobject
122 (
123 "collisionRecordsPairOrigProcOfOther",
125 ),
126 readOnProc
127 );
128 c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
129
130 labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
131 (
132 c.newIOobject
133 (
134 "collisionRecordsPairOrigIdOfOther",
136 ),
137 readOnProc
138 );
139 c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
140
141 pairDataFieldCompactIOField collisionRecordsPairData
142 (
143 c.newIOobject("collisionRecordsPairData", IOobject::MUST_READ),
144 readOnProc
145 );
146 c.checkFieldFieldIOobject(c, collisionRecordsPairData);
147
148 labelFieldCompactIOField collisionRecordsWallAccessed
149 (
150 c.newIOobject("collisionRecordsWallAccessed", IOobject::MUST_READ),
151 readOnProc
152 );
153 c.checkFieldFieldIOobject(c, collisionRecordsWallAccessed);
154
155 vectorFieldCompactIOField collisionRecordsWallPRel
156 (
157 c.newIOobject("collisionRecordsWallPRel", IOobject::MUST_READ),
158 readOnProc
159 );
160 c.checkFieldFieldIOobject(c, collisionRecordsWallPRel);
161
162 wallDataFieldCompactIOField collisionRecordsWallData
163 (
164 c.newIOobject("collisionRecordsWallData", IOobject::MUST_READ),
165 readOnProc
166 );
167 c.checkFieldFieldIOobject(c, collisionRecordsWallData);
168
169 label i = 0;
170
172 {
173 p.f_ = f[i];
174 p.angularMomentum_ = angularMomentum[i];
175 p.torque_ = torque[i];
176
177 p.collisionRecords_ = collisionRecordList
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
193template<class ParcelType>
194template<class CloudType>
196{
197 ParcelType::writeFields(c);
198
199 const label np = c.size();
200 const bool writeOnProc = c.size();
201
202
203 IOField<vector> f(c.newIOobject("f", IOobject::NO_READ), np);
204 IOField<vector> angMom
205 (
206 c.newIOobject("angularMomentum", IOobject::NO_READ),
207 np
208 );
209 IOField<vector> torque(c.newIOobject("torque", IOobject::NO_READ), np);
210
211 labelFieldCompactIOField collisionRecordsPairAccessed
212 (
213 c.newIOobject("collisionRecordsPairAccessed", IOobject::NO_READ),
214 np
215 );
216 labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
217 (
218 c.newIOobject
219 (
220 "collisionRecordsPairOrigProcOfOther",
222 ),
223 np
224 );
225 labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
226 (
227 c.newIOobject("collisionRecordsPairOrigIdOfOther", IOobject::NO_READ),
228 np
229 );
230 pairDataFieldCompactIOField collisionRecordsPairData
231 (
232 c.newIOobject("collisionRecordsPairData", IOobject::NO_READ),
233 np
234 );
235 labelFieldCompactIOField collisionRecordsWallAccessed
236 (
237 c.newIOobject("collisionRecordsWallAccessed", IOobject::NO_READ),
238 np
239 );
240 vectorFieldCompactIOField collisionRecordsWallPRel
241 (
242 c.newIOobject("collisionRecordsWallPRel", IOobject::NO_READ),
243 np
244 );
245 wallDataFieldCompactIOField collisionRecordsWallData
246 (
247 c.newIOobject("collisionRecordsWallData", IOobject::NO_READ),
248 np
249 );
250
251 label i = 0;
252 for (const CollidingParcel<ParcelType>& p : c)
253 {
254 f[i] = p.f();
255 angMom[i] = p.angularMomentum();
256 torque[i] = p.torque();
257
258 collisionRecordsPairAccessed[i] = p.collisionRecords().pairAccessed();
259 collisionRecordsPairOrigProcOfOther[i] =
260 p.collisionRecords().pairOrigProcOfOther();
261 collisionRecordsPairOrigIdOfOther[i] =
262 p.collisionRecords().pairOrigIdOfOther();
263 collisionRecordsPairData[i] = p.collisionRecords().pairData();
264 collisionRecordsWallAccessed[i] = p.collisionRecords().wallAccessed();
265 collisionRecordsWallPRel[i] = p.collisionRecords().wallPRel();
266 collisionRecordsWallData[i] = p.collisionRecords().wallData();
267
268 ++i;
269 }
270
271 f.write(writeOnProc);
272 angMom.write(writeOnProc);
273 torque.write(writeOnProc);
274
275 collisionRecordsPairAccessed.write(writeOnProc);
276 collisionRecordsPairOrigProcOfOther.write(writeOnProc);
277 collisionRecordsPairOrigIdOfOther.write(writeOnProc);
278 collisionRecordsPairData.write(writeOnProc);
279 collisionRecordsWallAccessed.write(writeOnProc);
280 collisionRecordsWallPRel.write(writeOnProc);
281 collisionRecordsWallData.write(writeOnProc);
282}
283
284
285template<class ParcelType>
287(
288 Ostream& os,
289 const wordRes& filters,
290 const word& delim,
291 const bool namesOnly
292) const
293{
294 ParcelType::writeProperties(os, filters, delim, namesOnly);
295
296 #undef writeProp
297 #define writeProp(Name, Value) \
298 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
299
300 writeProp("f", f_);
301 writeProp("angularMomentum", angularMomentum_);
302 writeProp("torque", torque_);
303 //writeProp("collisionRecords", collisionRecords_);
305 #undef writeProp
306}
307
308
309template<class ParcelType>
310template<class CloudType>
312(
313 CloudType& c,
314 const objectRegistry& obr
315)
316{
317 ParcelType::readObjects(c, obr);
318
319 if (!c.size()) return;
320
321 const auto& f = cloud::lookupIOField<vector>("f", obr);
322 const auto& angMom = cloud::lookupIOField<vector>("angularMomentum", obr);
323 const auto& torque = cloud::lookupIOField<vector>("torque", obr);
324
325 label i = 0;
327 {
328 p.f_ = f[i];
329 p.angularMomentum_ = angMom[i];
330 p.torque_ = torque[i];
331
332 ++i;
333 }
334}
335
336
337template<class ParcelType>
338template<class CloudType>
340(
341 const CloudType& c,
342 objectRegistry& obr
343)
344{
345 ParcelType::writeObjects(c, obr);
346
347 const label np = c.size();
348
349 auto& f = cloud::createIOField<vector>("f", np, obr);
350 auto& angMom = cloud::createIOField<vector>("angularMomentum", np, obr);
351 auto& torque = cloud::createIOField<vector>("torque", np, obr);
352
353 label i = 0;
354 for (const CollidingParcel<ParcelType>& p : c)
355 {
356 f[i] = p.f();
357 angMom[i] = p.angularMomentum();
358 torque[i] = p.torque();
359
360 ++i;
361 }
362}
363
364
365// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
366
367template<class ParcelType>
368Foam::Ostream& Foam::operator<<
369(
370 Ostream& os,
372)
373{
375 {
377 << token::SPACE << p.f_
378 << token::SPACE << p.angularMomentum_
379 << token::SPACE << p.torque_
380 << token::SPACE << p.collisionRecords_;
381 }
382 else
383 {
385 os.write
386 (
387 reinterpret_cast<const char*>(&p.f_),
389 );
390 os << p.collisionRecords();
391 }
392
394 return os;
395}
396
397
398// ************************************************************************* //
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 void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
const vector & angularMomentum() const
Return const access to angular momentum.
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.
const vector & f() const
Return const access to force.
collisionRecordList collisionRecords_
Particle collision records.
static void writeFields(const CloudType &c)
Write.
const vector & torque() const
Return const access to torque.
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.
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
streamFormat format() const noexcept
Get the current stream format.
@ ASCII
"ascii" (normal default)
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition IOstream.C:45
bool fatalCheckNativeSizes(const char *operation) const
Assert that the label/scalar byte-size associated with the stream are the native label/scalar sizes.
Definition IOstream.C:67
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
virtual Istream & read(token &)=0
Return next token from stream.
virtual Ostream & write(const char c) override
Write character.
Definition OBJstream.C:69
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static const IOField< Type > & lookupIOField(const word &fieldName, const objectRegistry &obr)
Lookup an IOField within object registry.
Definition cloud.H:202
static IOField< Type > & createIOField(const word &fieldName, const label nParticle, objectRegistry &obr)
Helper to construct IOField on a supplied object registry.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A class for handling character strings derived from std::string.
Definition string.H:76
@ SPACE
Space [isspace].
Definition token.H:144
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define FUNCTION_NAME
DSMCCloud< dsmcParcel > CloudType
vectorFieldCompactIOField pairDataFieldCompactIOField
CollisionRecordList< vector, vector > collisionRecordList
CompactIOField< labelField > labelFieldCompactIOField
Compact IO for a Field of labelField.
CompactIOField< vectorField > vectorFieldCompactIOField
Compact IO for a Field of vectorField.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vectorFieldCompactIOField wallDataFieldCompactIOField
labelList f(nPoints)