particleIO.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) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "particle.H"
30 #include "IOstreams.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
35 
36 const std::size_t Foam::particle::sizeofPosition
37 (
38  offsetof(particle, facei_) - offsetof(particle, coordinates_)
39 );
40 
41 const std::size_t Foam::particle::sizeofFields
42 (
43  sizeof(particle) - offsetof(particle, coordinates_)
44 );
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const polyMesh& mesh,
52  Istream& is,
53  const bool readFields,
54  const bool newFormat,
55  const bool doLocate
56 )
57 :
58  mesh_(mesh),
59  coordinates_(),
60  celli_(-1),
61  tetFacei_(-1),
62  tetPti_(-1),
63  facei_(-1),
64  stepFraction_(0.0),
65  origProc_(Pstream::myProcNo()),
66  origId_(-1)
67 {
69  readData(is, position, readFields, newFormat, doLocate);
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 (
77  Istream& is,
78  point& position,
79  const bool readFields,
80  const bool newFormat,
81  const bool doLocate
82 )
83 {
84  if (newFormat)
85  {
86  if (is.format() == IOstreamOption::ASCII)
87  {
88  is >> coordinates_ >> celli_ >> tetFacei_ >> tetPti_;
89  if (readFields)
90  {
91  is >> facei_ >> stepFraction_ >> origProc_ >> origId_;
92  }
93  }
94  else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
95  {
96  // Non-native label or scalar size
97 
98  is.beginRawRead();
99 
100  readRawScalar(is, coordinates_.data(), barycentric::nComponents);
101  readRawLabel(is, &celli_);
102  readRawLabel(is, &tetFacei_);
103  readRawLabel(is, &tetPti_);
104 
105  if (readFields)
106  {
107  readRawLabel(is, &facei_);
108  readRawScalar(is, &stepFraction_);
109  readRawLabel(is, &origProc_);
110  readRawLabel(is, &origId_);
111  }
112 
113  is.endRawRead();
114  }
115  else
116  {
117  if (readFields)
118  {
119  is.read(reinterpret_cast<char*>(&coordinates_), sizeofFields);
120  }
121  else
122  {
123  is.read(reinterpret_cast<char*>(&coordinates_), sizeofPosition);
124  }
125  }
126  }
127  else
128  {
129  positionsCompat1706 p;
130 
131  if (is.format() == IOstreamOption::ASCII)
132  {
133  is >> p.position >> p.celli;
134 
135  if (readFields)
136  {
137  is >> p.facei
138  >> p.stepFraction
139  >> p.tetFacei
140  >> p.tetPti
141  >> p.origProc
142  >> p.origId;
143  }
144  }
145  else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
146  {
147  // Non-native label or scalar size
148 
149  is.beginRawRead();
150 
151  readRawScalar(is, p.position.data(), vector::nComponents);
152  readRawLabel(is, &p.celli);
153 
154  if (readFields)
155  {
156  readRawLabel(is, &p.facei);
157  readRawScalar(is, &p.stepFraction);
158  readRawLabel(is, &p.tetFacei);
159  readRawLabel(is, &p.tetPti);
160  readRawLabel(is, &p.origProc);
161  readRawLabel(is, &p.origId);
162  }
163 
164  is.endRawRead();
165  }
166  else
167  {
168  if (readFields)
169  {
170  // Read whole struct
171  const size_t s =
172  (
173  sizeof(positionsCompat1706)
174  - offsetof(positionsCompat1706, position)
175  );
176  is.read(reinterpret_cast<char*>(&p.position), s);
177  }
178  else
179  {
180  // Read only position and cell
181  const size_t s =
182  (
183  offsetof(positionsCompat1706, facei)
184  - offsetof(positionsCompat1706, position)
185  );
186  is.read(reinterpret_cast<char*>(&p.position), s);
187  }
188  }
189 
190  if (readFields)
191  {
192  // Note: other position-based properties are set using locate(...)
193  stepFraction_ = p.stepFraction;
194  origProc_ = p.origProc;
195  origId_ = p.origId;
196  }
197 
198  // Preserve read position
199  position = p.position;
200 
201  if (doLocate)
202  {
203  locate
204  (
205  p.position,
206  nullptr,
207  p.celli,
208  false,
209  "Particle initialised with a location outside of the mesh."
210  );
211  }
212  }
214  // Check state of Istream
215  is.check(FUNCTION_NAME);
216 }
217 
218 
220 (
221  Ostream& os,
222  const wordRes& filters,
223  const word& delim,
224  const bool namesOnly
225 ) const
226 {
227  #undef writeProp
228  #define writeProp(Name, Value) \
229  particle::writeProperty(os, Name, Value, namesOnly, delim, filters)
230 
231  writeProp("coordinates", coordinates_);
232  writeProp("position", position());
233  writeProp("celli", celli_);
234  writeProp("tetFacei", tetFacei_);
235  writeProp("tetPti", tetPti_);
236  writeProp("facei", facei_);
237  writeProp("stepFraction", stepFraction_);
238  writeProp("origProc", origProc_);
239  writeProp("origId", origId_);
240 
241  #undef writeProp
242 }
243 
244 
246 {
248  {
249  os << coordinates_
250  << token::SPACE << celli_
251  << token::SPACE << tetFacei_
252  << token::SPACE << tetPti_;
253  }
254  else
255  {
256  os.write(reinterpret_cast<const char*>(&coordinates_), sizeofPosition);
257  }
258 
259  // Check state of Ostream
261 }
262 
263 
264 void Foam::particle::writePosition(Ostream& os) const
265 {
267  {
268  os << position() << token::SPACE << celli_;
269  }
270  else
271  {
272  positionsCompat1706 p;
273 
274  const size_t s =
275  (
276  offsetof(positionsCompat1706, facei)
277  - offsetof(positionsCompat1706, position)
278  );
279 
280  p.position = position();
281  p.celli = celli_;
282 
283  os.write(reinterpret_cast<const char*>(&p.position), s);
284  }
285 
286  // Check state of Ostream
288 }
289 
290 
291 Foam::Ostream& Foam::operator<<(Ostream& os, const particle& p)
292 {
294  {
295  os << p.coordinates_
296  << token::SPACE << p.celli_
297  << token::SPACE << p.tetFacei_
298  << token::SPACE << p.tetPti_
299  << token::SPACE << p.facei_
300  << token::SPACE << p.stepFraction_
301  << token::SPACE << p.origProc_
302  << token::SPACE << p.origId_;
303  }
304  else
305  {
306  os.write
307  (
308  reinterpret_cast<const char*>(&p.coordinates_),
309  particle::sizeofFields
310  );
311  }
312 
313  return os;
314 }
315 
316 
317 // ************************************************************************* //
virtual void writePosition(Ostream &os) const
Write the particle position and cell id.
Definition: particleIO.C:257
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const noexcept
Check if the scalar byte-size associated with the stream is the same as the given type...
Definition: IOstream.H:379
virtual Ostream & write(const char c) override
Write character.
Definition: OBJstream.C:69
static string propertyList_
String representation of properties.
Definition: particle.H:455
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
"ascii" (normal default)
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:272
virtual bool endRawRead()=0
End of low-level raw binary read.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Base particle class.
Definition: particle.H:69
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
void readData(Istream &is, point &position, const bool readFields, const bool newFormat, const bool doLocate)
Read particle from stream. Optionally (for old format) return.
Definition: particleIO.C:69
virtual Istream & read(token &)=0
Return next token from stream.
dynamicFvMesh & mesh
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
A class for handling words, derived from Foam::string.
Definition: word.H:63
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition: particle.C:507
Inter-processor communications stream.
Definition: Pstream.H:57
Space [isspace].
Definition: token.H:131
void writeCoordinates(Ostream &os) const
Write the particle barycentric coordinates and cell info.
Definition: particleIO.C:238
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:109
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
label readRawLabel(Istream &is)
Read raw label from binary stream.
Definition: label.C:39
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
virtual bool beginRawRead()=0
Start of low-level raw binary read.
#define writeProp(Name, Value)
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const noexcept
Check if the label byte-size associated with the stream is the same as the given type.
Definition: IOstream.H:368
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual particle properties to stream.
Definition: particleIO.C:213
A class for handling character strings derived from std::string.
Definition: string.H:72
streamFormat format() const noexcept
Get the current stream format.
vector position() const
Return current particle position.
Definition: particleI.H:283
static string propertyList()
Definition: particle.H:455