moleculeIO.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2019 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 "molecule.H"
30#include "IOstreams.H"
31#include "moleculeCloud.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35const std::size_t Foam::molecule::sizeofFields
36(
37 offsetof(molecule, siteForces_) - offsetof(molecule, Q_)
38);
39
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
44(
45 const polyMesh& mesh,
46 Istream& is,
47 bool readFields,
48 bool newFormat
49)
50:
51 particle(mesh, is, readFields, newFormat),
52 Q_(Zero),
53 v_(Zero),
54 a_(Zero),
55 pi_(Zero),
56 tau_(Zero),
57 specialPosition_(Zero),
58 potentialEnergy_(0.0),
59 rf_(Zero),
60 special_(0),
61 id_(0),
62 siteForces_(),
63 sitePositions_()
64{
65 if (readFields)
66 {
67 if (is.format() == IOstream::ASCII)
68 {
69 is >> Q_
70 >> v_
71 >> a_
72 >> pi_
73 >> tau_
74 >> specialPosition_
75 >> potentialEnergy_
76 >> rf_
77 >> special_
78 >> id_;
79 }
80 else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
81 {
82 // Non-native label or scalar size
83
84 is.beginRawRead();
85
86 readRawScalar(is, Q_.data(), tensor::nComponents);
87 readRawScalar(is, v_.data(), vector::nComponents);
88 readRawScalar(is, a_.data(), vector::nComponents);
89 readRawScalar(is, pi_.data(), vector::nComponents);
90 readRawScalar(is, tau_.data(), vector::nComponents);
91 readRawScalar(is, specialPosition_.data(), vector::nComponents);
92 readRawScalar(is, &potentialEnergy_);
93 readRawScalar(is, rf_.data(), tensor::nComponents);
94 readRawLabel(is, &special_);
95 readRawLabel(is, &id_);
96
97 is.endRawRead();
98 }
99 else
100 {
101 is.read(reinterpret_cast<char*>(&Q_), sizeofFields);
102 }
103
104 is >> siteForces_ >> sitePositions_;
105 }
106
108}
109
110
112{
113 const bool valid = mC.size();
114
116
118 mC.checkFieldIOobject(mC, Q);
119
121 mC.checkFieldIOobject(mC, v);
122
124 mC.checkFieldIOobject(mC, a);
125
127 mC.checkFieldIOobject(mC, pi);
128
130 mC.checkFieldIOobject(mC, tau);
131
132 IOField<vector> specialPosition
133 (
134 mC.fieldIOobject("specialPosition", IOobject::MUST_READ),
135 valid
136 );
137 mC.checkFieldIOobject(mC, specialPosition);
138
139 IOField<label> special
140 (
141 mC.fieldIOobject("special", IOobject::MUST_READ),
142 valid
143 );
144 mC.checkFieldIOobject(mC, special);
145
147 mC.checkFieldIOobject(mC, id);
148
149 label i = 0;
150 for (molecule& mol : mC)
151 {
152 mol.Q_ = Q[i];
153 mol.v_ = v[i];
154 mol.a_ = a[i];
155 mol.pi_ = pi[i];
156 mol.tau_ = tau[i];
157 mol.specialPosition_ = specialPosition[i];
158 mol.special_ = special[i];
159 mol.id_ = id[i];
160
161 ++i;
162 }
163}
164
165
167{
169
170 const label np = mC.size();
171 const bool valid = np;
172
178 IOField<vector> specialPosition
179 (
180 mC.fieldIOobject("specialPosition", IOobject::NO_READ),
181 np
182 );
183 IOField<label> special(mC.fieldIOobject("special", IOobject::NO_READ), np);
185
186 // Post processing fields
187
188 IOField<vector> piGlobal
189 (
190 mC.fieldIOobject("piGlobal", IOobject::NO_READ),
191 np
192 );
193
194 IOField<vector> tauGlobal
195 (
196 mC.fieldIOobject("tauGlobal", IOobject::NO_READ),
197 np
198 );
199
200 IOField<vector> orientation1
201 (
202 mC.fieldIOobject("orientation1", IOobject::NO_READ),
203 np
204 );
205
206 IOField<vector> orientation2
207 (
208 mC.fieldIOobject("orientation2", IOobject::NO_READ),
209 np
210 );
211
212 IOField<vector> orientation3
213 (
214 mC.fieldIOobject("orientation3", IOobject::NO_READ),
215 np
216 );
217
218 label i = 0;
219 for (const molecule& mol : mC)
220 {
221 Q[i] = mol.Q_;
222 v[i] = mol.v_;
223 a[i] = mol.a_;
224 pi[i] = mol.pi_;
225 tau[i] = mol.tau_;
226 specialPosition[i] = mol.specialPosition_;
227 special[i] = mol.special_;
228 id[i] = mol.id_;
229
230 piGlobal[i] = mol.Q_ & mol.pi_;
231 tauGlobal[i] = mol.Q_ & mol.tau_;
232
233 orientation1[i] = mol.Q_ & vector(1,0,0);
234 orientation2[i] = mol.Q_ & vector(0,1,0);
235 orientation3[i] = mol.Q_ & vector(0,0,1);
236
237 ++i;
238 }
239
240 Q.write(valid);
241 v.write(valid);
242 a.write(valid);
243 pi.write(valid);
244 tau.write(valid);
245 specialPosition.write(valid);
246 special.write(valid);
247 id.write(valid);
248
249 piGlobal.write(valid);
250 tauGlobal.write(valid);
251
252 orientation1.write(valid);
253 orientation2.write(valid);
254 orientation3.write(valid);
255
256 Info<< "writeFields " << mC.name() << endl;
257
258 if (isA<moleculeCloud>(mC))
259 {
260 const moleculeCloud& m = dynamic_cast<const moleculeCloud&>(mC);
261
262 m.writeXYZ
263 (
264 m.mesh().time().timePath()/cloud::prefix/"moleculeCloud.xmol"
265 );
266 }
267}
268
269
270// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
271
273{
274 if (os.format() == IOstream::ASCII)
275 {
276 os << token::SPACE << static_cast<const particle&>(mol)
277 << token::SPACE << mol.Q_
278 << token::SPACE << mol.v_
279 << token::SPACE << mol.a_
280 << token::SPACE << mol.pi_
281 << token::SPACE << mol.tau_
282 << token::SPACE << mol.specialPosition_
283 << token::SPACE << mol.potentialEnergy_
284 << token::SPACE << mol.rf_
285 << token::SPACE << mol.special_
286 << token::SPACE << mol.id_
287 << token::SPACE << mol.siteForces_
288 << token::SPACE << mol.sitePositions_;
289 }
290 else
291 {
292 os << static_cast<const particle&>(mol);
293 os.write
294 (
295 reinterpret_cast<const char*>(&mol.Q_),
297 );
298 os << mol.siteForces_ << mol.sitePositions_;
299 }
300
302 return os;
303}
304
305
306// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Base cloud calls templated on particle type.
Definition: Cloud.H:68
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:211
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:191
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
streamFormat format() const noexcept
Get the current stream format.
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const noexcept
Definition: IOstream.H:300
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const noexcept
Definition: IOstream.H:309
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
virtual bool endRawRead()=0
End of low-level raw binary read.
virtual bool beginRawRead()=0
Start of low-level raw binary read.
virtual Istream & read(token &)=0
Return next token from stream.
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
fileName timePath() const
Return current time path.
Definition: Time.H:375
Cmpt * data() noexcept
Return pointer to the first data element.
Definition: VectorSpaceI.H:185
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:87
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
const polyMesh & mesh() const
void writeXYZ(const fileName &fName) const
Write molecule sites in XYZ format.
Foam::molecule.
Definition: molecule.H:70
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: molecule.H:74
static void readFields(Cloud< molecule > &mC)
Definition: moleculeIO.C:111
const Time & time() const noexcept
Return time registry.
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
Base particle class.
Definition: particle.H:79
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual bool write(const bool valid=true) const
Write using setting from DB.
@ SPACE
Space [isspace].
Definition: token.H:125
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label readRawLabel(Istream &is)
Read raw label from binary stream.
Definition: label.C:46
Vector< scalar > vector
Definition: vector.H:61