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.newIOobject("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> UCorrect(c.newIOobject("UCorrect", IOobject::NO_READ), np);
122 
123  label i = 0;
124 
125  for (const MPPICParcel<ParcelType>& p : c)
126  {
127  UCorrect[i] = p.UCorrect();
128 
129  ++i;
130  }
132  UCorrect.write(writeOnProc);
133 }
134 
135 
136 template<class ParcelType>
138 (
139  Ostream& os,
140  const wordRes& filters,
141  const word& delim,
142  const bool namesOnly
143 ) const
144 {
145  ParcelType::writeProperties(os, filters, delim, namesOnly);
146 
147  #undef writeProp
148  #define writeProp(Name, Value) \
149  ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
150 
151  writeProp("UCorrect", UCorrect_);
152 
153  #undef writeProp
154 }
155 
156 
157 template<class ParcelType>
158 template<class CloudType>
160 (
161  CloudType& c,
162  const objectRegistry& obr
163 )
164 {
165  ParcelType::readObjects(c, obr);
166 
167  // const bool readOnProc = c.size();
168  if (!c.size()) return;
169 
170  const auto& UCorrect = cloud::lookupIOField<vector>("UCorrect", obr);
171 
172  label i = 0;
173  for (MPPICParcel<ParcelType>& p : c)
174  {
175  p.UCorrect() = UCorrect[i];
176 
177  ++i;
178  }
179 }
180 
181 
182 template<class ParcelType>
183 template<class CloudType>
185 (
186  const CloudType& c,
187  objectRegistry& obr
188 )
189 {
190  ParcelType::writeObjects(c, obr);
191 
192  const label np = c.size();
193  // const bool writeOnProc = c.size();
194 
195  auto& UCorrect = cloud::createIOField<vector>("UCorrect", np, obr);
196 
197  label i = 0;
198  for (const MPPICParcel<ParcelType>& p : c)
199  {
200  UCorrect[i] = p.UCorrect();
201 
202  ++i;
203  }
204 }
205 
206 
207 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
208 
209 template<class ParcelType>
210 Foam::Ostream& Foam::operator<<
211 (
212  Ostream& os,
213  const MPPICParcel<ParcelType>& p
214 )
215 {
217  {
218  os << static_cast<const ParcelType&>(p)
219  << token::SPACE << p.UCorrect();
220  }
221  else
222  {
223  os << static_cast<const ParcelType&>(p);
224  os.write
225  (
226  reinterpret_cast<const char*>(&p.UCorrect_),
228  );
229  }
230 
232  return os;
233 }
234 
235 
236 // ************************************************************************* //
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:75
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