MPPICParcelIO.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) 2013-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 "MPPICParcel.H"
30 #include "IOstreams.H"
31 #include "IOField.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 template<class ParcelType>
38 
39 
40 template<class ParcelType>
42 (
43  sizeof(MPPICParcel<ParcelType>) - sizeof(ParcelType)
44 );
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 template<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  UCorrect_(Zero)
60 {
61  if (readFields)
62  {
63  if (is.format() == IOstreamOption::ASCII)
64  {
65  is >> UCorrect_;
66  }
67  else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
68  {
69  // Non-native label or scalar size
70 
71  is.beginRawRead();
72 
73  readRawScalar(is, UCorrect_.data(), vector::nComponents);
74 
75  is.endRawRead();
76  }
77  else
78  {
79  is.read(reinterpret_cast<char*>(&UCorrect_), sizeofFields);
80  }
81  }
82 
83  is.check(FUNCTION_NAME);
84 }
85 
86 
87 template<class ParcelType>
88 template<class CloudType>
90 {
91  const bool readOnProc = c.size();
92 
94 
95  IOField<vector> UCorrect
96  (
97  c.fieldIOobject("UCorrect", IOobject::MUST_READ),
98  readOnProc
99  );
100  c.checkFieldIOobject(c, UCorrect);
101 
102  label i = 0;
103  for (MPPICParcel<ParcelType>& p : c)
104  {
105  p.UCorrect_ = UCorrect[i];
106 
107  ++i;
108  }
109 }
110 
111 
112 template<class ParcelType>
113 template<class CloudType>
115 {
117 
118  const label np = c.size();
119  const bool writeOnProc = c.size();
120 
121  IOField<vector>
122  UCorrect(c.fieldIOobject("UCorrect", IOobject::NO_READ), np);
123 
124  label i = 0;
125 
126  for (const MPPICParcel<ParcelType>& p : c)
127  {
128  UCorrect[i] = p.UCorrect();
129 
130  ++i;
131  }
133  UCorrect.write(writeOnProc);
134 }
135 
136 
137 template<class ParcelType>
139 (
140  Ostream& os,
141  const wordRes& filters,
142  const word& delim,
143  const bool namesOnly
144 ) const
145 {
146  ParcelType::writeProperties(os, filters, delim, namesOnly);
147 
148  #undef writeProp
149  #define writeProp(Name, Value) \
150  ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
151 
152  writeProp("UCorrect", UCorrect_);
153 
154  #undef writeProp
155 }
156 
157 
158 template<class ParcelType>
159 template<class CloudType>
161 (
162  CloudType& c,
163  const objectRegistry& obr
164 )
165 {
166  ParcelType::readObjects(c, obr);
167 
168  // const bool readOnProc = c.size();
169  if (!c.size()) return;
170 
171  const auto& UCorrect = cloud::lookupIOField<vector>("UCorrect", obr);
172 
173  label i = 0;
174  for (MPPICParcel<ParcelType>& p : c)
175  {
176  p.UCorrect() = UCorrect[i];
177 
178  ++i;
179  }
180 }
181 
182 
183 template<class ParcelType>
184 template<class CloudType>
186 (
187  const CloudType& c,
188  objectRegistry& obr
189 )
190 {
191  ParcelType::writeObjects(c, obr);
192 
193  const label np = c.size();
194  // const bool writeOnProc = c.size();
195 
196  auto& UCorrect = cloud::createIOField<vector>("UCorrect", np, obr);
197 
198  label i = 0;
199  for (const MPPICParcel<ParcelType>& p : c)
200  {
201  UCorrect[i] = p.UCorrect();
202 
203  ++i;
204  }
205 }
206 
207 
208 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
209 
210 template<class ParcelType>
211 Foam::Ostream& Foam::operator<<
212 (
213  Ostream& os,
214  const MPPICParcel<ParcelType>& p
215 )
216 {
218  {
219  os << static_cast<const ParcelType&>(p)
220  << token::SPACE << p.UCorrect();
221  }
222  else
223  {
224  os << static_cast<const ParcelType&>(p);
225  os.write
226  (
227  reinterpret_cast<const char*>(&p.UCorrect_),
229  );
230  }
231 
233  return os;
234 }
235 
236 
237 // ************************************************************************* //
DSMCCloud< dsmcParcel > CloudType
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
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)
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
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
#define writeProp(Name, Value)
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
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
Wrapper around kinematic parcel types to add MPPIC modelling.
Definition: MPPICParcel.H:54
Space [isspace].
Definition: token.H:131
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
static void writeFields(const CloudType &c)
Write.
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: MPPICParcel.H:82
vector UCorrect_
Velocity correction due to collisions [m/s].
Definition: MPPICParcel.H:196
MPPICParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: MPPICParcelI.H:25
virtual bool beginRawRead()=0
Start of low-level raw binary read.
const dimensionedScalar c
Speed of light in a vacuum.
Nothing to be read.
static void readFields(CloudType &c)
Read.
Definition: MPPICParcelIO.C:82
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
Registry of regIOobjects.
A class for handling character strings derived from std::string.
Definition: string.H:72
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
A primitive field of type <T> with automated input and output.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
streamFormat format() const noexcept
Get the current stream format.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127