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 | Copyright (C) 2013-2016 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  if (!c.size())
81  {
82  return;
83  }
84 
86 
87  IOField<vector> UCorrect(c.fieldIOobject("UCorrect", IOobject::MUST_READ));
88  c.checkFieldIOobject(c, UCorrect);
89 
90  label i = 0;
91 
92  forAllIter(typename CloudType, c, iter)
93  {
94  MPPICParcel<ParcelType>& p = iter();
95 
96  p.UCorrect_ = UCorrect[i];
97 
98  i++;
99  }
100 }
101 
102 
103 template<class ParcelType>
104 template<class CloudType>
106 {
107  ParcelType::writeFields(c);
108 
109  label np = c.size();
110 
112  UCorrect(c.fieldIOobject("UCorrect", IOobject::NO_READ), np);
113 
114  label i = 0;
115 
116  forAllConstIter(typename CloudType, c, iter)
117  {
118  const MPPICParcel<ParcelType>& p = iter();
119 
120  UCorrect[i] = p.UCorrect();
121 
122  i++;
123  }
124 
125  UCorrect.write();
126 }
127 
128 
129 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
130 
131 template<class ParcelType>
132 Foam::Ostream& Foam::operator<<
133 (
134  Ostream& os,
136 )
137 {
138  if (os.format() == IOstream::ASCII)
139  {
140  os << static_cast<const ParcelType&>(p)
141  << token::SPACE << p.UCorrect();
142  }
143  else
144  {
145  os << static_cast<const ParcelType&>(p);
146  os.write
147  (
148  reinterpret_cast<const char*>(&p.UCorrect_),
150  );
151  }
152 
153  os.check
154  (
155  "Ostream& operator<<(Ostream&, const MPPICParcel<ParcelType>&)"
156  );
157 
158  return os;
159 }
160 
161 
162 // ************************************************************************* //
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:215
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
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
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 kinematic parcel types to add MPPIC modelling.
Definition: MPPICParcel.H:52
static const zero Zero
Definition: zero.H:91
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static void writeFields(const CloudType &c)
Write.
label size() const
Definition: Cloud.H:175
MPPICParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
Definition: MPPICParcelI.H:30
vector UCorrect_
Velocity correction due to collisions [m/s].
Definition: MPPICParcel.H:163
const dimensionedScalar c
Speed of light in a vacuum.
const vector & UCorrect() const
Return const access to correction velocity.
Definition: MPPICParcelI.H:81
virtual bool write() const
Write using setting from DB.
static void readFields(CloudType &c)
Read.
Definition: MPPICParcelIO.C:78
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:195
volScalarField & p
A class for handling character strings derived from std::string.
Definition: string.H:74
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68