ThermoParcelIO.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) 2011-2018 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 "ThermoParcel.H"
27 #include "IOstreams.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<class ParcelType>
32 Foam::string Foam::ThermoParcel<ParcelType>::propertyList_ =
34 
35 template<class ParcelType>
37 (
38  sizeof(ThermoParcel<ParcelType>)
39  - offsetof(ThermoParcel<ParcelType>, T_)
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  T_(0.0),
55  Cp_(0.0)
56 {
57  if (readFields)
58  {
59  if (is.format() == IOstream::ASCII)
60  {
61  T_ = readScalar(is);
62  Cp_ = readScalar(is);
63  }
64  else
65  {
66  is.read(reinterpret_cast<char*>(&T_), sizeofFields_);
67  }
68  }
69 
70  // Check state of Istream
71  is.check
72  (
73  "ThermoParcel::ThermoParcel(const polyMesh&, Istream&, bool)"
74  );
75 }
76 
77 
78 template<class ParcelType>
79 template<class CloudType>
81 {
82  bool valid = c.size();
83 
85 
86  IOField<scalar> T(c.fieldIOobject("T", IOobject::MUST_READ), valid);
87  c.checkFieldIOobject(c, T);
88 
89  IOField<scalar> Cp(c.fieldIOobject("Cp", IOobject::MUST_READ), valid);
90  c.checkFieldIOobject(c, Cp);
91 
92 
93  label i = 0;
94  forAllIter(typename Cloud<ThermoParcel<ParcelType>>, c, iter)
95  {
96  ThermoParcel<ParcelType>& p = iter();
97 
98  p.T_ = T[i];
99  p.Cp_ = Cp[i];
100  i++;
101  }
102 }
103 
104 
105 template<class ParcelType>
106 template<class CloudType>
108 {
109  ParcelType::writeFields(c);
110 
111  label np = c.size();
112 
113  IOField<scalar> T(c.fieldIOobject("T", IOobject::NO_READ), np);
114  IOField<scalar> Cp(c.fieldIOobject("Cp", IOobject::NO_READ), np);
115 
116  label i = 0;
118  {
119  const ThermoParcel<ParcelType>& p = iter();
120 
121  T[i] = p.T_;
122  Cp[i] = p.Cp_;
123  i++;
124  }
125 
126  T.write(np > 0);
127  Cp.write(np > 0);
128 }
129 
130 
131 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
132 
133 template<class ParcelType>
134 Foam::Ostream& Foam::operator<<
135 (
136  Ostream& os,
138 )
139 {
140  if (os.format() == IOstream::ASCII)
141  {
142  os << static_cast<const ParcelType&>(p)
143  << token::SPACE << p.T()
144  << token::SPACE << p.Cp();
145  }
146  else
147  {
148  os << static_cast<const ParcelType&>(p);
149  os.write
150  (
151  reinterpret_cast<const char*>(&p.T_),
153  );
154  }
155 
156  // Check state of Ostream
157  os.check
158  (
159  "Ostream& operator<<(Ostream&, const ThermoParcel<ParcelType>&)"
160  );
161 
162  return os;
163 }
164 
165 
166 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
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 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
scalar Cp_
Specific heat capacity [J/kg/K].
Definition: ThermoParcel.H:252
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
#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
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: ThermoParcelI.H:75
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:151
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
scalar T_
Temperature [K].
Definition: ThermoParcel.H:249
virtual Istream & read(token &)=0
Return next token from stream.
static void readFields(CloudType &c)
Read.
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Base cloud calls templated on particle type.
Definition: Cloud.H:52
static void writeFields(const CloudType &c)
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const volScalarField & T
const dimensionedScalar c
Speed of light in a vacuum.
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:187
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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:69
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:167
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
Definition: ThermoParcel.H:51