All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
particleIO.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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, coordinates_)
36 );
37 
38 const std::size_t Foam::particle::sizeofFields_
39 (
40  sizeof(particle) - offsetof(particle, coordinates_)
41 );
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 :
48  mesh_(mesh),
49  coordinates_(),
50  celli_(-1),
51  tetFacei_(-1),
52  tetPti_(-1),
53  facei_(-1),
54  stepFraction_(0.0),
55  behind_(0.0),
56  nBehind_(0),
57  origProc_(Pstream::myProcNo()),
58  origId_(-1)
59 {
60  if (is.format() == IOstream::ASCII)
61  {
62  is >> coordinates_ >> celli_ >> tetFacei_ >> tetPti_;
63 
64  if (readFields)
65  {
66  is >> facei_ >> stepFraction_ >> behind_ >> nBehind_
67  >> origProc_ >> origId_;
68  }
69  }
70  else
71  {
72  if (readFields)
73  {
74  is.read(reinterpret_cast<char*>(&coordinates_), sizeofFields_);
75  }
76  else
77  {
78  is.read(reinterpret_cast<char*>(&coordinates_), sizeofPosition_);
79  }
80  }
81 
82  // Check state of Istream
83  is.check("particle::particle(Istream&, bool)");
84 }
85 
86 
88 {
89  if (os.format() == IOstream::ASCII)
90  {
91  os << coordinates_
92  << token::SPACE << celli_
93  << token::SPACE << tetFacei_
94  << token::SPACE << tetPti_;
95  }
96  else
97  {
98  os.write(reinterpret_cast<const char*>(&coordinates_), sizeofPosition_);
99  }
100 
101  // Check state of Ostream
102  os.check("particle::writePosition(Ostream& os, bool) const");
103 }
104 
105 
107 {
108  if (os.format() == IOstream::ASCII)
109  {
110  os << p.coordinates_
111  << token::SPACE << p.celli_
112  << token::SPACE << p.tetFacei_
113  << token::SPACE << p.tetPti_
114  << token::SPACE << p.facei_
115  << token::SPACE << p.stepFraction_
116  << token::SPACE << p.behind_
117  << token::SPACE << p.nBehind_
118  << token::SPACE << p.origProc_
119  << token::SPACE << p.origId_;
120  }
121  else
122  {
123  os.write
124  (
125  reinterpret_cast<const char*>(&p.coordinates_),
126  particle::sizeofFields_
127  );
128  }
129 
130  return os;
131 }
132 
133 
134 // ************************************************************************* //
virtual Ostream & write(const char)=0
Write character.
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
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Base particle class.
Definition: particle.H:84
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
virtual Istream & read(token &)=0
Return next token from stream.
dynamicFvMesh & mesh
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition: particle.C:512
Inter-processor communications stream.
Definition: Pstream.H:53
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writePosition(Ostream &) const
Write the particle position and cell.
Definition: particleIO.C:87
Ostream & operator<<(Ostream &, const ensightPart &)
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:76
static string propertyList()
Definition: particle.H:368