Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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 "MPPICParcel.H"
30#include "IOstreams.H"
31#include "IOField.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35template<class ParcelType>
38
39
40template<class ParcelType>
42(
43 sizeof(MPPICParcel<ParcelType>) - sizeof(ParcelType)
44);
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49template<class ParcelType>
51(
52 const polyMesh& mesh,
53 Istream& is,
54 bool readFields,
55 bool newFormat
56)
57:
58 ParcelType(mesh, is, readFields, newFormat),
60{
61 if (readFields)
62 {
63 if (is.format() == IOstreamOption::ASCII)
64 {
65 is >> UCorrect_;
66 }
67 else
68 {
69 // No non-native streaming
71
72 is.read(reinterpret_cast<char*>(&UCorrect_), sizeofFields);
73 }
74 }
77}
78
79
80template<class ParcelType>
81template<class CloudType>
83{
84 const bool readOnProc = c.size();
85
86 ParcelType::readFields(c);
87
88 IOField<vector> UCorrect
89 (
90 c.newIOobject("UCorrect", IOobject::MUST_READ),
91 readOnProc
92 );
93 c.checkFieldIOobject(c, UCorrect);
94
95 label i = 0;
96 for (MPPICParcel<ParcelType>& p : c)
97 {
98 p.UCorrect_ = UCorrect[i];
99
100 ++i;
101 }
102}
103
104
105template<class ParcelType>
106template<class CloudType>
108{
109 ParcelType::writeFields(c);
110
111 const label np = c.size();
112 const bool writeOnProc = c.size();
113
114 IOField<vector> UCorrect(c.newIOobject("UCorrect", IOobject::NO_READ), np);
115
116 label i = 0;
117
118 for (const MPPICParcel<ParcelType>& p : c)
119 {
120 UCorrect[i] = p.UCorrect();
121
122 ++i;
124
125 UCorrect.write(writeOnProc);
126}
127
128
129template<class ParcelType>
131(
132 Ostream& os,
133 const wordRes& filters,
134 const word& delim,
135 const bool namesOnly
136) const
137{
138 ParcelType::writeProperties(os, filters, delim, namesOnly);
139
140 #undef writeProp
141 #define writeProp(Name, Value) \
142 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
143
144 writeProp("UCorrect", UCorrect_);
146 #undef writeProp
147}
148
149
150template<class ParcelType>
151template<class CloudType>
153(
154 CloudType& c,
155 const objectRegistry& obr
156)
157{
158 ParcelType::readObjects(c, obr);
159
160 // const bool readOnProc = c.size();
161 if (!c.size()) return;
162
163 const auto& UCorrect = cloud::lookupIOField<vector>("UCorrect", obr);
164
165 label i = 0;
166 for (MPPICParcel<ParcelType>& p : c)
167 {
168 p.UCorrect() = UCorrect[i];
169
170 ++i;
171 }
172}
173
174
175template<class ParcelType>
176template<class CloudType>
178(
179 const CloudType& c,
180 objectRegistry& obr
181)
182{
183 ParcelType::writeObjects(c, obr);
184
185 const label np = c.size();
186 // const bool writeOnProc = c.size();
187
188 auto& UCorrect = cloud::createIOField<vector>("UCorrect", np, obr);
189
190 label i = 0;
191 for (const MPPICParcel<ParcelType>& p : c)
192 {
193 UCorrect[i] = p.UCorrect();
194
195 ++i;
196 }
197}
198
199
200// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
201
202template<class ParcelType>
203Foam::Ostream& Foam::operator<<
204(
205 Ostream& os,
207)
208{
210 {
212 << token::SPACE << p.UCorrect();
213 }
214 else
215 {
217 os.write
218 (
219 reinterpret_cast<const char*>(&p.UCorrect_),
221 );
222 }
223
225 return os;
226}
227
228
229// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
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.
Wrapper around kinematic parcel types to add MPPIC modelling.
Definition MPPICParcel.H:76
const vector & UCorrect() const
Return const access to correction velocity.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition MPPICParcel.H:82
MPPICParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
static void writeFields(const CloudType &c)
Write.
vector UCorrect_
Velocity correction due to collisions [m/s].
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
static void readFields(CloudType &c)
Read.
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
static const IOField< Type > & lookupIOField(const word &fieldName, const objectRegistry &obr)
Lookup an IOField within object registry.
Definition cloud.H:202
static IOField< Type > & createIOField(const word &fieldName, const label nParticle, objectRegistry &obr)
Helper to construct IOField on a supplied object registry.
Registry of regIOobjects.
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
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define FUNCTION_NAME
DSMCCloud< dsmcParcel > CloudType
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127