Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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 "particle.H"
30#include "IOstreams.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
35
36const std::size_t Foam::particle::sizeofPosition
37(
38 offsetof(particle, facei_) - offsetof(particle, coordinates_)
39);
40
41const std::size_t Foam::particle::sizeofFields
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
95 {
96 // No non-native streaming
98
99 if (readFields)
100 {
101 is.read(reinterpret_cast<char*>(&coordinates_), sizeofFields);
102 }
103 else
104 {
105 is.read(reinterpret_cast<char*>(&coordinates_), sizeofPosition);
106 }
107 }
108 }
109 else
110 {
111 positionsCompat1706 p;
112
113 if (is.format() == IOstreamOption::ASCII)
114 {
115 is >> p.position >> p.celli;
116
117 if (readFields)
118 {
119 is >> p.facei
120 >> p.stepFraction
121 >> p.tetFacei
122 >> p.tetPti
123 >> p.origProc
124 >> p.origId;
125 }
126 }
127 else
128 {
129 // No non-native streaming
131
132 if (readFields)
133 {
134 // Read whole struct
135 const size_t s =
136 (
137 sizeof(positionsCompat1706)
138 - offsetof(positionsCompat1706, position)
139 );
140 is.read(reinterpret_cast<char*>(&p.position), s);
141 }
142 else
143 {
144 // Read only position and cell
145 const size_t s =
146 (
147 offsetof(positionsCompat1706, facei)
148 - offsetof(positionsCompat1706, position)
149 );
150 is.read(reinterpret_cast<char*>(&p.position), s);
151 }
152 }
153
154 if (readFields)
155 {
156 // Note: other position-based properties are set using locate(...)
157 stepFraction_ = p.stepFraction;
158 origProc_ = p.origProc;
159 origId_ = p.origId;
160 }
161
162 // Preserve read position
163 position = p.position;
164
165 if (doLocate)
166 {
167 locate
168 (
169 p.position,
170 nullptr,
171 p.celli,
172 false,
173 "Particle initialised with a location outside of the mesh."
174 );
175 }
177
178 // Check state of Istream
180}
181
182
184(
185 Ostream& os,
186 const wordRes& filters,
187 const word& delim,
188 const bool namesOnly
189) const
190{
191 #undef writeProp
192 #define writeProp(Name, Value) \
193 particle::writeProperty(os, Name, Value, namesOnly, delim, filters)
194
195 writeProp("coordinates", coordinates_);
196 writeProp("position", position());
197 writeProp("celli", celli_);
198 writeProp("tetFacei", tetFacei_);
199 writeProp("tetPti", tetPti_);
200 writeProp("facei", facei_);
201 writeProp("stepFraction", stepFraction_);
202 writeProp("origProc", origProc_);
203 writeProp("origId", origId_);
204
205 #undef writeProp
206}
207
208
210{
211 if (os.format() == IOstreamOption::ASCII)
212 {
213 os << coordinates_
214 << token::SPACE << celli_
215 << token::SPACE << tetFacei_
216 << token::SPACE << tetPti_;
217 }
218 else
219 {
220 os.write(reinterpret_cast<const char*>(&coordinates_), sizeofPosition);
222
223 // Check state of Ostream
224 os.check(FUNCTION_NAME);
225}
226
227
228void Foam::particle::writePosition(Ostream& os) const
229{
230 if (os.format() == IOstreamOption::ASCII)
231 {
232 os << position() << token::SPACE << celli_;
233 }
234 else
235 {
236 positionsCompat1706 p;
237
238 const size_t s =
239 (
240 offsetof(positionsCompat1706, facei)
241 - offsetof(positionsCompat1706, position)
242 );
243
244 p.position = position();
245 p.celli = celli_;
246
247 os.write(reinterpret_cast<const char*>(&p.position), s);
249
250 // Check state of Ostream
251 os.check(FUNCTION_NAME);
252}
253
254
255Foam::Ostream& Foam::operator<<(Ostream& os, const particle& p)
256{
258 {
259 os << p.coordinates_
260 << token::SPACE << p.celli_
261 << token::SPACE << p.tetFacei_
262 << token::SPACE << p.tetPti_
263 << token::SPACE << p.facei_
264 << token::SPACE << p.stepFraction_
265 << token::SPACE << p.origProc_
266 << token::SPACE << p.origId_;
267 }
268 else
269 {
270 os.write
271 (
272 reinterpret_cast<const char*>(&p.coordinates_),
273 particle::sizeofFields
274 );
275 }
276
277 return os;
278}
279
280
281// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
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
Inter-processor communications stream.
Definition Pstream.H:59
Base particle class.
Definition particle.H:72
vector position() const
Return current particle position.
Definition particleI.H:283
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
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
static string propertyList()
Definition particle.H:455
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual particle properties to stream.
Definition particleIO.C:177
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition particle.C:507
static string propertyList_
String representation of properties.
Definition particle.H:455
virtual void writePosition(Ostream &os) const
Write the particle position and cell id.
Definition particleIO.C:221
void writeCoordinates(Ostream &os) const
Write the particle barycentric coordinates and cell info.
Definition particleIO.C:202
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
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
OBJstream os(runTime.globalPath()/outputName)
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))
#define writeProp(Name, Value)
#define FUNCTION_NAME
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.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
vector point
Point is a vector.
Definition point.H:37
Old particle positions content for OpenFOAM-1706 and earlier.
Definition particle.H:123