All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
MPPICParcelIO.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) 2013-2020 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 "MPPICParcel.H"
27 #include "IOstreams.H"
28 #include "IOField.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 template<class ParcelType>
33 Foam::string Foam::MPPICParcel<ParcelType>::propertyList_ =
35 
36 template<class ParcelType>
38 (
39  sizeof(MPPICParcel<ParcelType>) - sizeof(ParcelType)
40 );
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 template<class ParcelType>
47 (
48  const polyMesh& mesh,
49  Istream& is,
50  bool readFields
51 )
52 :
53  ParcelType(mesh, is, readFields),
54  UCorrect_(Zero)
55 {
56  if (readFields)
57  {
58  if (is.format() == IOstream::ASCII)
59  {
60  is >> UCorrect_;
61  }
62  else
63  {
64  is.read(reinterpret_cast<char*>(&UCorrect_), sizeofFields_);
65  }
66  }
67 
68  is.check
69  (
70  "MPPICParcel<ParcelType>::Collisions"
71  "(const polyMesh&, Istream&, bool)"
72  );
73 }
74 
75 
76 template<class ParcelType>
77 template<class CloudType>
79 {
80  bool valid = c.size();
81 
83 
84  IOField<vector> UCorrect
85  (
86  c.fieldIOobject("UCorrect", IOobject::MUST_READ),
87  valid
88  );
89  c.checkFieldIOobject(c, UCorrect);
90 
91  label i = 0;
92 
93  forAllIter(typename CloudType, c, iter)
94  {
95  MPPICParcel<ParcelType>& p = iter();
96 
97  p.UCorrect_ = UCorrect[i];
98 
99  i++;
100  }
101 }
102 
103 
104 template<class ParcelType>
105 template<class CloudType>
107 {
108  ParcelType::writeFields(c);
109 
110  label np = c.size();
111 
113  UCorrect(c.fieldIOobject("UCorrect", IOobject::NO_READ), np);
114 
115  label i = 0;
116 
117  forAllConstIter(typename CloudType, c, iter)
118  {
119  const MPPICParcel<ParcelType>& p = iter();
120 
121  UCorrect[i] = p.UCorrect();
122 
123  i++;
124  }
125 
126  UCorrect.write(np > 0);
127 }
128 
129 
130 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
131 
132 template<class ParcelType>
133 Foam::Ostream& Foam::operator<<
134 (
135  Ostream& os,
137 )
138 {
139  if (os.format() == IOstream::ASCII)
140  {
141  os << static_cast<const ParcelType&>(p)
142  << token::SPACE << p.UCorrect();
143  }
144  else
145  {
146  os << static_cast<const ParcelType&>(p);
147  os.write
148  (
149  reinterpret_cast<const char*>(&p.UCorrect_),
151  );
152  }
153 
154  os.check
155  (
156  "Ostream& operator<<(Ostream&, const MPPICParcel<ParcelType>&)"
157  );
158 
159  return os;
160 }
161 
162 
163 // ************************************************************************* //
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
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:190
const dimensionedScalar c
Speed of light in a vacuum.
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.
Wrapper around parcel types to add MPPIC modelling.
Definition: MPPICParcel.H:52
static const zero Zero
Definition: zero.H:97
streamFormat format() const
Return current stream format.
Definition: IOstream.H:374
const vector & UCorrect() const
Return const access to correction velocity.
Definition: MPPICParcelI.H:59
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static void writeFields(const CloudType &c)
Write.
vector UCorrect_
Velocity correction due to collisions [m/s].
Definition: MPPICParcel.H:168
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:30
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:190
static void readFields(CloudType &c)
Read.
Definition: MPPICParcelIO.C:78
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
virtual bool write(const bool write=true) const
Write using setting from DB.
volScalarField & p
A class for handling character strings derived from std::string.
Definition: string.H:76
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:170
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50