ReactingParcel.H
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 Class
25  Foam::ReactingParcel
26 
27 Description
28  Reacting parcel class with one/two-way coupling with the continuous
29  phase.
30 
31 SourceFiles
32  ReactingParcelI.H
33  ReactingParcel.C
34  ReactingParcelIO.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef ReactingParcel_H
39 #define ReactingParcel_H
40 
41 #include "particle.H"
42 #include "interpolation.H"
43 #include "fluidThermo.H"
44 #include "demandDrivenEntry.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 template<class ParcelType>
52 class ReactingParcel;
53 
54 template<class ParcelType>
55 Ostream& operator<<
56 (
57  Ostream&,
59 );
60 
61 
62 /*---------------------------------------------------------------------------*\
63  Class ReactingParcelName Declaration
64 \*---------------------------------------------------------------------------*/
65 
67 
68 
69 /*---------------------------------------------------------------------------*\
70  Class ReactingParcel Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 template<class ParcelType>
74 class ReactingParcel
75 :
76  public ParcelType,
77  public ReactingParcelName
78 {
79  // Private Data
80 
81  //- Size in bytes of the fields
82  static const std::size_t sizeofFields_;
83 
84 
85 public:
86 
87  //- Class to hold reacting parcel constant properties
88  class constantProperties
89  :
90  public ParcelType::constantProperties
91  {
92  // Private Data
93 
94  //- Minimum pressure [Pa]
96 
97  //- Constant volume flag - e.g. during mass transfer
98  demandDrivenEntry<bool> constantVolume_;
99 
100 
101  public:
102 
103  // Constructors
104 
105  //- Null constructor
107 
108  //- Copy constructor
110 
111  //- Construct from dictionary
112  constantProperties(const dictionary& parentDict);
113 
114 
115  // Access
116 
117  //- Return const access to the minimum pressure
118  inline scalar pMin() const;
119 
120  //- Return const access to the constant volume flag
121  inline bool constantVolume() const;
122  };
123 
124 
125  //- Use base tracking data
126  typedef typename ParcelType::trackingData trackingData;
127 
128 
129 protected:
130 
131  // Protected data
132 
133  // Parcel properties
134 
135  //- Mass fractions of mixture []
136  scalarField Y_;
137 
138 
139  // Protected Member Functions
140 
141  //- Calculate Phase change
142  template<class TrackCloudType>
143  void calcPhaseChange
144  (
145  TrackCloudType& cloud,
146  trackingData& td,
147  const scalar dt, // timestep
148  const scalar Re, // Reynolds number
149  const scalar Pr, // Prandtl number
150  const scalar Ts, // Surface temperature
151  const scalar nus, // Surface kinematic viscosity
152  const scalar d, // diameter
153  const scalar T, // temperature
154  const scalar mass, // mass
155  const label idPhase, // id of phase involved in phase change
156  const scalar YPhase, // total mass fraction
157  const scalarField& YComponents, // component mass fractions
158  scalarField& dMassPC, // mass transfer - local to parcel
159  scalar& Sh, // explicit parcel enthalpy source
160  scalar& N, // flux of species emitted from parcel
161  scalar& NCpW, // sum of N*Cp*W of emission species
162  scalarField& Cs // carrier conc. of emission species
163  );
164 
165  //- Update mass fraction
166  scalar updateMassFraction
167  (
168  const scalar mass0,
169  const scalarField& dMass,
170  scalarField& Y
171  ) const;
172 
173 
174 public:
175 
176  // Static Data Members
177 
178  //- String representation of properties
180  (
181  ParcelType,
182  " nPhases(Y1..YN)"
183  );
184 
185 
186  // Constructors
187 
188  //- Construct from mesh, coordinates and topology
189  // Other properties initialised as null
190  inline ReactingParcel
191  (
192  const polyMesh& mesh,
193  const barycentric& coordinates,
194  const label celli,
195  const label tetFacei,
196  const label tetPti,
197  const label facei
198  );
199 
200  //- Construct from a position and a cell, searching for the rest of the
201  // required topology. Other properties are initialised as null.
202  inline ReactingParcel
203  (
204  const polyMesh& mesh,
205  const vector& position,
206  const label celli
207  );
208 
209  //- Construct from Istream
210  ReactingParcel(Istream& is, bool readFields = true);
211 
212  //- Construct as a copy
214 
215  //- Construct and return a clone
216  virtual autoPtr<particle> clone() const
217  {
219  }
220 
221  //- Construct from Istream and return
223  {
224  return autoPtr<ReactingParcel>(new ReactingParcel(is));
225  }
226 
227 
228  // Member Functions
229 
230  // Access
231 
232  //- Return const access to mass fractions of mixture []
233  inline const scalarField& Y() const;
234 
235 
236  // Edit
237 
238  //- Return access to mass fractions of mixture []
239  inline scalarField& Y();
240 
241 
242  // Main calculation loop
243 
244  //- Set cell values
245  template<class TrackCloudType>
246  void setCellValues(TrackCloudType& cloud, trackingData& td);
247 
248  //- Correct cell values using latest transfer information
249  template<class TrackCloudType>
251  (
252  TrackCloudType& cloud,
253  trackingData& td,
254  const scalar dt
255  );
256 
257  //- Correct surface values due to emitted species
258  template<class TrackCloudType>
260  (
261  TrackCloudType& cloud,
262  trackingData& td,
263  const scalar T,
264  const scalarField& Cs,
265  scalar& rhos,
266  scalar& mus,
267  scalar& Prs,
268  scalar& kappas
269  );
270 
271  //- Update parcel properties over the time interval
272  template<class TrackCloudType>
273  void calc
274  (
275  TrackCloudType& cloud,
276  trackingData& td,
277  const scalar dt
278  );
279 
280 
281  // I-O
282 
283  //- Read
284  template<class CloudType, class CompositionType>
285  static void readFields
286  (
287  CloudType& c,
288  const CompositionType& compModel
289  );
290 
291  //- Read - no composition
292  template<class CloudType>
293  static void readFields(CloudType& c);
294 
295  //- Write
296  template<class CloudType, class CompositionType>
297  static void writeFields
298  (
299  const CloudType& c,
300  const CompositionType& compModel
301  );
302 
303  //- Write - composition supplied
304  template<class CloudType>
305  static void writeFields(const CloudType& c);
306 
307 
308  // Ostream Operator
309 
310  friend Ostream& operator<< <ParcelType>
311  (
312  Ostream&,
314  );
315 };
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 } // End namespace Foam
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 #include "ReactingParcelI.H"
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 #ifdef NoRepository
329  #include "ReactingParcel.C"
330 #endif
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
334 #endif
335 
336 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Class to hold reacting parcel constant properties.
bool constantVolume() const
Return const access to the constant volume flag.
scalar pMin() const
Return const access to the minimum pressure.
Reacting parcel class with one/two-way coupling with the continuous phase.
static autoPtr< ReactingParcel > New(Istream &is)
Construct from Istream and return.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
void calcPhaseChange(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar Ts, const scalar nus, const scalar d, const scalar T, const scalar mass, const label idPhase, const scalar YPhase, const scalarField &YComponents, scalarField &dMassPC, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs)
Calculate Phase change.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
const scalarField & Y() const
Return const access to mass fractions of mixture [].
void correctSurfaceValues(TrackCloudType &cloud, trackingData &td, const scalar T, const scalarField &Cs, scalar &rhos, scalar &mus, scalar &Prs, scalar &kappas)
Correct surface values due to emitted species.
scalarField Y_
Mass fractions of mixture [].
ParcelType::trackingData trackingData
Use base tracking data.
scalar updateMassFraction(const scalar mass0, const scalarField &dMass, scalarField &Y) const
Update mass fraction.
virtual autoPtr< particle > clone() const
Construct and return a clone.
AddToPropertyList(ParcelType, " nPhases(Y1..YN)")
String representation of properties.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
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.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
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
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:753
TemplateName(FvFaceCellWave)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
volScalarField & p