ReactingParcelIO.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-2023 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 "ReactingParcel.H"
27 #include "IOstreams.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<class ParcelType>
34 
35 template<class ParcelType>
37 (
38  0
39 );
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class ParcelType>
46 :
47  ParcelType(is, readFields),
48  Y_(0)
49 {
50  if (readFields)
51  {
53 
54  is >> Ymix;
55 
56  Y_.transfer(Ymix);
57  }
58 
59  // Check state of Istream
60  is.check
61  (
62  "ReactingParcel<ParcelType>::ReactingParcel"
63  "("
64  "const polyMesh&, "
65  "Istream&, "
66  "bool"
67  ")"
68  );
69 }
70 
71 
72 template<class ParcelType>
73 template<class CloudType>
75 {
77 }
78 
79 
80 template<class ParcelType>
81 template<class CloudType, class CompositionType>
83 (
84  CloudType& c,
85  const CompositionType& compModel
86 )
87 {
88  bool valid = c.size();
89 
91 
92  // Get names and sizes for each Y...
93  const wordList& phaseTypes = compModel.phaseTypes();
94  const label nPhases = phaseTypes.size();
95  wordList stateLabels(nPhases, "");
96  if (compModel.nPhase() == 1)
97  {
98  stateLabels = compModel.stateLabels()[0];
99  }
100 
101  // Set storage for each Y... for each parcel
102  forAllIter(typename CloudType, c, iter)
103  {
104  ReactingParcel<ParcelType>& p = iter();
105  p.Y_.setSize(nPhases, 0.0);
106  }
107 
108  // Populate Y for each parcel
109  forAll(phaseTypes, j)
110  {
112  (
113  c.fieldIOobject
114  (
115  "Y" + phaseTypes[j] + stateLabels[j],
117  ),
118  valid
119  );
120 
121  label i = 0;
122  forAllIter(typename CloudType, c, iter)
123  {
124  ReactingParcel<ParcelType>& p = iter();
125  p.Y_[j] = Y[i++];
126  }
127  }
128 }
129 
130 
131 template<class ParcelType>
132 template<class CloudType>
134 {
135  ParcelType::writeFields(c);
136 }
137 
138 
139 template<class ParcelType>
140 template<class CloudType, class CompositionType>
142 (
143  const CloudType& c,
144  const CompositionType& compModel
145 )
146 {
147  ParcelType::writeFields(c, compModel);
148 
149  const label np = c.size();
150 
151  {
152  // Write the composition fractions
153  const wordList& phaseTypes = compModel.phaseTypes();
154  wordList stateLabels(phaseTypes.size(), "");
155  if (compModel.nPhase() == 1)
156  {
157  stateLabels = compModel.stateLabels()[0];
158  }
159 
160  forAll(phaseTypes, j)
161  {
163  (
164  c.fieldIOobject
165  (
166  "Y" + phaseTypes[j] + stateLabels[j],
168  ),
169  np
170  );
171  label i = 0;
172  forAllConstIter(typename CloudType, c, iter)
173  {
174  const ReactingParcel<ParcelType>& p = iter();
175  Y[i++] = p.Y()[j];
176  }
177 
178  Y.write(np > 0);
179  }
180  }
181 }
182 
183 
184 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
185 
186 template<class ParcelType>
187 Foam::Ostream& Foam::operator<<
188 (
189  Ostream& os,
191 )
192 {
193  if (os.format() == IOstream::ASCII)
194  {
195  os << static_cast<const ParcelType&>(p)
196  << token::SPACE << p.Y();
197  }
198  else
199  {
200  os << static_cast<const ParcelType&>(p);
201  os << p.Y();
202  }
203 
204  // Check state of Ostream
205  os.check
206  (
207  "Ostream& operator<<(Ostream&, const ReactingParcel<ParcelType>&)"
208  );
209 
210  return os;
211 }
212 
213 
214 // ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Reacting parcel class with one/two-way coupling with the continuous phase.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
scalarField Y_
Mass fractions of mixture [].
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
A class for handling character strings derived from std::string.
Definition: string.H:79
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar c
Speed of light in a vacuum.
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
volScalarField & p
PtrList< volScalarField > & Y