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 | Copyright (C) 2011-2016 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 "particle.H"
27 #include "IOstreams.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 Foam::string Foam::particle::propertyList_ = Foam::particle::propertyList();
32 
33 const std::size_t Foam::particle::sizeofPosition_
34 (
35  offsetof(particle, facei_) - offsetof(particle, position_)
36 );
37 
38 const std::size_t Foam::particle::sizeofFields_
39 (
40  sizeof(particle) - offsetof(particle, position_)
41 );
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 :
48  mesh_(mesh),
49  position_(),
50  celli_(-1),
51  facei_(-1),
52  stepFraction_(0.0),
53  tetFacei_(-1),
54  tetPtI_(-1),
55  origProc_(Pstream::myProcNo()),
56  origId_(-1)
57 {
58  if (is.format() == IOstream::ASCII)
59  {
60  is >> position_ >> celli_;
61 
62  if (readFields)
63  {
64  is >> facei_
65  >> stepFraction_
66  >> tetFacei_
67  >> tetPtI_
68  >> origProc_
69  >> origId_;
70  }
71  }
72  else
73  {
74  if (readFields)
75  {
76  is.read(reinterpret_cast<char*>(&position_), sizeofFields_);
77  }
78  else
79  {
80  is.read(reinterpret_cast<char*>(&position_), sizeofPosition_);
81  }
82  }
83 
84  // Check state of Istream
85  is.check("particle::particle(Istream&, bool)");
86 }
87 
88 
90 {
91  if (os.format() == IOstream::ASCII)
92  {
93  os << position_ << token::SPACE << celli_;
94  }
95  else
96  {
97  os.write(reinterpret_cast<const char*>(&position_), sizeofPosition_);
98  }
99 
100  // Check state of Ostream
101  os.check("particle::writePosition(Ostream& os, bool) const");
102 }
103 
104 
106 {
107  if (os.format() == IOstream::ASCII)
108  {
109  os << p.position_
110  << token::SPACE << p.celli_
111  << token::SPACE << p.facei_
113  << token::SPACE << p.tetFacei_
114  << token::SPACE << p.tetPtI_
115  << token::SPACE << p.origProc_
116  << token::SPACE << p.origId_;
117  }
118  else
119  {
120  os.write
121  (
122  reinterpret_cast<const char*>(&p.position_),
123  particle::sizeofFields_
124  );
125  }
126 
127  return os;
128 }
129 
130 
131 // ************************************************************************* //
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label tetPtI_
Index of the point on the face that defines the decomposed.
Definition: particle.H:158
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
particle(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from components.
Definition: particle.C:48
vector position_
Position of particle.
Definition: particle.H:140
label origId_
Local particle id on originating processor.
Definition: particle.H:164
Base particle class.
Definition: particle.H:78
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
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
label origProc_
Originating processor id.
Definition: particle.H:161
virtual Istream & read(token &)=0
Return next token from stream.
dynamicFvMesh & mesh
Inter-processor communications stream.
Definition: Pstream.H:53
scalar stepFraction_
Fraction of time-step completed.
Definition: particle.H:149
label celli_
Index of the cell it is in.
Definition: particle.H:143
label tetFacei_
Index of the face that owns the decomposed tet that the.
Definition: particle.H:153
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void writePosition(Ostream &) const
Write the particle position and cell.
Definition: particleIO.C:89
Ostream & operator<<(Ostream &, const ensightPart &)
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
A class for handling character strings derived from std::string.
Definition: string.H:74
static string propertyList()
Definition: particle.H:319
label facei_
Face index if the particle is on a face otherwise -1.
Definition: particle.H:146